Intrepid
http://trilinos.sandia.gov/packages/docs/r10.10/packages/intrepid/test/Discretization/Basis/HGRAD_WEDGE_C2_FEM/test_01.cpp
00001 // @HEADER
00002 // ************************************************************************
00003 //
00004 //                           Intrepid Package
00005 //                 Copyright (2007) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Pavel Bochev  (pbboche@sandia.gov)
00038 //                    Denis Ridzal  (dridzal@sandia.gov), or
00039 //                    Kara Peterson (kjpeter@sandia.gov)
00040 //
00041 // ************************************************************************
00042 // @HEADER
00043 
00048 #include "Intrepid_FieldContainer.hpp"
00049 #include "Intrepid_HGRAD_WEDGE_C2_FEM.hpp"
00050 #include "Teuchos_oblackholestream.hpp"
00051 #include "Teuchos_RCP.hpp"
00052 #include "Teuchos_GlobalMPISession.hpp"
00053 
00054 using namespace std;
00055 using namespace Intrepid;
00056 
00057 #define INTREPID_TEST_COMMAND( S , throwCounter, nException )                                                              \
00058 {                                                                                                                          \
00059   ++nException;                                                                                                            \
00060   try {                                                                                                                    \
00061     S ;                                                                                                                    \
00062   }                                                                                                                        \
00063   catch (std::logic_error err) {                                                                                           \
00064       ++throwCounter;                                                                                                      \
00065       *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
00066       *outStream << err.what() << '\n';                                                                                    \
00067       *outStream << "-------------------------------------------------------------------------------" << "\n\n";           \
00068   };                                                                                                                       \
00069 }
00070 
00071 int main(int argc, char *argv[]) {
00072   
00073   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00074 
00075   // This little trick lets us print to std::cout only if
00076   // a (dummy) command-line argument is provided.
00077   int iprint     = argc - 1;
00078   Teuchos::RCP<std::ostream> outStream;
00079   Teuchos::oblackholestream bhs; // outputs nothing
00080   if (iprint > 0)
00081     outStream = Teuchos::rcp(&std::cout, false);
00082   else
00083     outStream = Teuchos::rcp(&bhs, false);
00084   
00085   // Save the format state of the original std::cout.
00086   Teuchos::oblackholestream oldFormatState;
00087   oldFormatState.copyfmt(std::cout);
00088   
00089   *outStream \
00090     << "===============================================================================\n" \
00091     << "|                                                                             |\n" \
00092     << "|                 Unit Test (Basis_HGRAD_WEDGE_C2_FEM)                        |\n" \
00093     << "|                                                                             |\n" \
00094     << "|     1) Conversion of Dof tags into Dof ordinals and back                    |\n" \
00095     << "|     2) Basis values for VALUE, GRAD, and Dk operators                       |\n" \
00096     << "|                                                                             |\n" \
00097     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00098     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00099     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00100     << "|                                                                             |\n" \
00101     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00102     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00103     << "|                                                                             |\n" \
00104     << "===============================================================================\n"\
00105     << "| TEST 1: Basis creation, exception testing                                   |\n"\
00106     << "===============================================================================\n";
00107   
00108   // Define basis and error flag
00109   Basis_HGRAD_WEDGE_C2_FEM<double, FieldContainer<double> > wedgeBasis;
00110   int errorFlag = 0;
00111 
00112   // Initialize throw counter for exception testing
00113   int nException     = 0;
00114   int throwCounter   = 0;
00115 
00116   // Nodes of Wedde<18>: vertices, edge midpoints, Quadrilateral face centers 
00117   FieldContainer<double> wedgeNodes(18, 3);
00118   wedgeNodes(0,0) =  0.0;  wedgeNodes(0,1) =  0.0;  wedgeNodes(0,2) = -1.0;  
00119   wedgeNodes(1,0) =  1.0;  wedgeNodes(1,1) =  0.0;  wedgeNodes(1,2) = -1.0;  
00120   wedgeNodes(2,0) =  0.0;  wedgeNodes(2,1) =  1.0;  wedgeNodes(2,2) = -1.0;
00121   wedgeNodes(3,0) =  0.0;  wedgeNodes(3,1) =  0.0;  wedgeNodes(3,2) =  1.0;  
00122   wedgeNodes(4,0) =  1.0;  wedgeNodes(4,1) =  0.0;  wedgeNodes(4,2) =  1.0;  
00123   wedgeNodes(5,0) =  0.0;  wedgeNodes(5,1) =  1.0;  wedgeNodes(5,2) =  1.0;
00124     
00125   wedgeNodes(6,0) =  0.5;  wedgeNodes(6,1) =  0.0;  wedgeNodes(6,2) = -1.0;  
00126   wedgeNodes(7,0) =  0.5;  wedgeNodes(7,1) =  0.5;  wedgeNodes(7,2) = -1.0;  
00127   wedgeNodes(8,0) =  0.0;  wedgeNodes(8,1) =  0.5;  wedgeNodes(8,2) = -1.0;
00128   wedgeNodes(9,0) =  0.0;  wedgeNodes(9,1) =  0.0;  wedgeNodes(9,2) =  0.0;
00129   wedgeNodes(10,0)=  1.0;  wedgeNodes(10,1)=  0.0;  wedgeNodes(10,2)=  0.0;  
00130   wedgeNodes(11,0)=  0.0;  wedgeNodes(11,1)=  1.0;  wedgeNodes(11,2)=  0.0;  
00131  
00132   wedgeNodes(12,0)=  0.5;  wedgeNodes(12,1)=  0.0;  wedgeNodes(12,2)=  1.0;  
00133   wedgeNodes(13,0)=  0.5;  wedgeNodes(13,1)=  0.5;  wedgeNodes(13,2)=  1.0;  
00134   wedgeNodes(14,0)=  0.0;  wedgeNodes(14,1)=  0.5;  wedgeNodes(14,2)=  1.0;  
00135   wedgeNodes(15,0)=  0.5;  wedgeNodes(15,1)=  0.0;  wedgeNodes(15,2)=  0.0;  
00136   wedgeNodes(16,0)=  0.5;  wedgeNodes(16,1)=  0.5;  wedgeNodes(16,2)=  0.0;  
00137   wedgeNodes(17,0)=  0.0;  wedgeNodes(17,1)=  0.5;  wedgeNodes(17,2)=  0.0;  
00138 
00139 
00140   // Generic array for the output values; needs to be properly resized depending on the operator type
00141   FieldContainer<double> vals;
00142 
00143   try{
00144     // exception #1: CURL cannot be applied to scalar functions
00145     // resize vals to rank-3 container with dimensions (num. points, num. basis functions)
00146     vals.resize(wedgeBasis.getCardinality(), wedgeNodes.dimension(0), 3 );
00147     INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_DIV), throwCounter, nException );
00148 
00149     // exception #2: DIV cannot be applied to scalar functions
00150     // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
00151     vals.resize(wedgeBasis.getCardinality(), wedgeNodes.dimension(0) );
00152     INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_DIV), throwCounter, nException );
00153         
00154     // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 
00155     // getDofTag() to access invalid array elements thereby causing bounds check exception
00156     // exception #3
00157     INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(3,0,0), throwCounter, nException );
00158     // exception #4
00159     INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(1,1,1), throwCounter, nException );
00160     // exception #5
00161     INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(0,9,0), throwCounter, nException );
00162     // exception #6
00163     INTREPID_TEST_COMMAND( wedgeBasis.getDofTag(18), throwCounter, nException );
00164     // exception #7
00165     INTREPID_TEST_COMMAND( wedgeBasis.getDofTag(-1), throwCounter, nException );
00166     
00167 #ifdef HAVE_INTREPID_DEBUG
00168     // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays
00169     // exception #8: input points array must be of rank-2
00170     FieldContainer<double> badPoints1(4, 5, 3);
00171     INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00172     
00173     // exception #9 dimension 1 in the input point array must equal space dimension of the cell
00174     FieldContainer<double> badPoints2(4, wedgeBasis.getBaseCellTopology().getDimension() + 1);
00175     INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00176     
00177     // exception #10 output values must be of rank-2 for OPERATOR_VALUE
00178     FieldContainer<double> badVals1(4, 3, 1);
00179     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals1, wedgeNodes, OPERATOR_VALUE), throwCounter, nException );
00180     
00181     // exception #11 output values must be of rank-3 for OPERATOR_GRAD
00182     FieldContainer<double> badVals2(4, 3);
00183     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals2, wedgeNodes, OPERATOR_GRAD), throwCounter, nException );
00184     
00185     // exception #12 output values must be of rank-3 for OPERATOR_D1
00186     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals2, wedgeNodes, OPERATOR_D1), throwCounter, nException );
00187     
00188     // exception #13 output values must be of rank-3 for OPERATOR_D2
00189     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals2, wedgeNodes, OPERATOR_D2), throwCounter, nException );
00190     
00191     // exception #14 incorrect 0th dimension of output array (must equal number of basis functions)
00192     FieldContainer<double> badVals3(wedgeBasis.getCardinality() + 1, wedgeNodes.dimension(0));
00193     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals3, wedgeNodes, OPERATOR_VALUE), throwCounter, nException );
00194     
00195     // exception #15 incorrect 1st dimension of output array (must equal number of points)
00196     FieldContainer<double> badVals4(wedgeBasis.getCardinality(), wedgeNodes.dimension(0) + 1);
00197     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals4, wedgeNodes, OPERATOR_VALUE), throwCounter, nException );
00198     
00199     // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
00200     FieldContainer<double> badVals5(wedgeBasis.getCardinality(), wedgeNodes.dimension(0), wedgeBasis.getBaseCellTopology().getDimension() - 1);
00201     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals5, wedgeNodes, OPERATOR_GRAD), throwCounter, nException );
00202     
00203     // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 3D)
00204     FieldContainer<double> badVals6(wedgeBasis.getCardinality(), wedgeNodes.dimension(0), 40);
00205     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals6, wedgeNodes, OPERATOR_D2), throwCounter, nException );
00206     
00207     // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 3D)
00208     INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals6, wedgeNodes, OPERATOR_D3), throwCounter, nException );
00209 #endif
00210 
00211   }
00212   catch (std::logic_error err) {
00213     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00214     *outStream << err.what() << '\n';
00215     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00216     errorFlag = -1000;
00217   };
00218   
00219   // Check if number of thrown exceptions matches the one we expect - 18
00220   if (throwCounter != nException) {
00221     errorFlag++;
00222     *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00223   }
00224   
00225   *outStream \
00226     << "\n"
00227     << "===============================================================================\n"\
00228     << "| TEST 2: correctness of tag to enum and enum to tag lookups                  |\n"\
00229     << "===============================================================================\n";
00230   
00231   try{
00232     std::vector<std::vector<int> > allTags = wedgeBasis.getAllDofTags();
00233     
00234     // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
00235     for (unsigned i = 0; i < allTags.size(); i++) {
00236       int bfOrd  = wedgeBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
00237       
00238       std::vector<int> myTag = wedgeBasis.getDofTag(bfOrd);
00239        if( !( (myTag[0] == allTags[i][0]) &&
00240               (myTag[1] == allTags[i][1]) &&
00241               (myTag[2] == allTags[i][2]) &&
00242               (myTag[3] == allTags[i][3]) ) ) {
00243         errorFlag++;
00244         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00245         *outStream << " getDofOrdinal( {" 
00246           << allTags[i][0] << ", " 
00247           << allTags[i][1] << ", " 
00248           << allTags[i][2] << ", " 
00249           << allTags[i][3] << "}) = " << bfOrd <<" but \n";   
00250         *outStream << " getDofTag(" << bfOrd << ") = { "
00251           << myTag[0] << ", " 
00252           << myTag[1] << ", " 
00253           << myTag[2] << ", " 
00254           << myTag[3] << "}\n";        
00255       }
00256     }
00257     
00258     // Now do the same but loop over basis functions
00259     for( int bfOrd = 0; bfOrd < wedgeBasis.getCardinality(); bfOrd++) {
00260       std::vector<int> myTag  = wedgeBasis.getDofTag(bfOrd);
00261       int myBfOrd = wedgeBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
00262       if( bfOrd != myBfOrd) {
00263         errorFlag++;
00264         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00265         *outStream << " getDofTag(" << bfOrd << ") = { "
00266           << myTag[0] << ", " 
00267           << myTag[1] << ", " 
00268           << myTag[2] << ", " 
00269           << myTag[3] << "} but getDofOrdinal({" 
00270           << myTag[0] << ", " 
00271           << myTag[1] << ", " 
00272           << myTag[2] << ", " 
00273           << myTag[3] << "} ) = " << myBfOrd << "\n";
00274       }
00275     }
00276   }
00277   catch (std::logic_error err){
00278     *outStream << err.what() << "\n\n";
00279     errorFlag = -1000;
00280   };
00281   
00282   *outStream \
00283     << "\n"
00284     << "===============================================================================\n"\
00285     << "| TEST 3: correctness of basis function values                                |\n"\
00286     << "===============================================================================\n";
00287   
00288   outStream -> precision(20);
00289   
00290   // VALUE: correct basis function values in (F,P) format
00291   double basisValues[] = {
00292     1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, \
00293     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, \
00294     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, \
00295     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00296     0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00297     0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00298     0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00299     1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, \
00300     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, \
00301     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, \
00302     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00303     0, 0, 0, 0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00304     0, 0, 0, 0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00305     0, 0, 1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00306     1.00, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00  };
00307   
00308   // GRAD, D1, D2, D3 and D4 test values are stored in files due to their large size
00309   std::string     fileName;
00310   std::ifstream   dataFile;
00311   
00312   // GRAD and D1 values are stored in (F,P,D) format in a data file. Read file and do the test
00313   std::vector<double> basisGrads;           // Flat array for the gradient values.
00314   
00315   fileName = "./testdata/WEDGE_C2_GradVals.dat";
00316   dataFile.open(fileName.c_str());
00317   TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00318                       ">>> ERROR (HGRAD_WEDGE_C2/test01): could not open GRAD values data file, test aborted.");
00319   while (!dataFile.eof() ){
00320     double temp;
00321     string line;                            // string for one line of input file
00322     std::getline(dataFile, line);           // get next line from file
00323     stringstream data_line(line);           // convert to stringstream
00324     while(data_line >> temp){               // extract value from line
00325       basisGrads.push_back(temp);           // push into vector
00326     }
00327   }
00328   // It turns out that just closing and then opening the ifstream variable does not reset it
00329   // and subsequent open() command fails. One fix is to explicitely clear the ifstream, or
00330   // scope the variables.
00331   dataFile.close();
00332   dataFile.clear();
00333   
00334   
00335   //D2: flat array with the values of D2 applied to basis functions. Multi-index is (F,P,D2cardinality)
00336   std::vector<double> basisD2;
00337   
00338   fileName = "./testdata/WEDGE_C2_D2Vals.dat";
00339   dataFile.open(fileName.c_str());
00340   TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00341                       ">>> ERROR (HGRAD_WEDGE_C2/test01): could not open D2 values data file, test aborted.");
00342   
00343   while (!dataFile.eof() ){
00344     double temp;
00345     string line;                            // string for one line of input file
00346     std::getline(dataFile, line);           // get next line from file
00347     stringstream data_line(line);           // convert to stringstream
00348     while(data_line >> temp){               // extract value from line
00349       basisD2.push_back(temp);              // push into vector
00350     }
00351   }
00352   dataFile.close();
00353   dataFile.clear();
00354   
00355   
00356   //D3: flat array with the values of D3 applied to basis functions. Multi-index is (F,P,D3cardinality)
00357   std::vector<double> basisD3;
00358   
00359   fileName = "./testdata/WEDGE_C2_D3Vals.dat";
00360   dataFile.open(fileName.c_str());
00361   TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00362                       ">>> ERROR (HGRAD_WEDGE_C2/test01): could not open D3 values data file, test aborted.");
00363   
00364   while (!dataFile.eof() ){
00365     double temp;
00366     string line;                            // string for one line of input file
00367     std::getline(dataFile, line);           // get next line from file
00368     stringstream data_line(line);           // convert to stringstream
00369     while(data_line >> temp){               // extract value from line
00370       basisD3.push_back(temp);              // push into vector
00371     }
00372   }
00373   dataFile.close();
00374   dataFile.clear();
00375   
00376   
00377   //D4: flat array with the values of D3 applied to basis functions. Multi-index is (F,P,D4cardinality)
00378   std::vector<double> basisD4;
00379   
00380   fileName = "./testdata/WEDGE_C2_D4Vals.dat";
00381   dataFile.open(fileName.c_str());
00382   TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00383                       ">>> ERROR (HGRAD_WEDGE_C2/test01): could not open D4 values data file, test aborted.");
00384   
00385   while (!dataFile.eof() ){
00386     double temp;
00387     string line;                            // string for one line of input file
00388     std::getline(dataFile, line);           // get next line from file
00389     stringstream data_line(line);           // convert to stringstream
00390     while(data_line >> temp){               // extract value from line
00391       basisD4.push_back(temp);              // push into vector
00392     }
00393   }
00394   dataFile.close();
00395   dataFile.clear();
00396 
00397   
00398   try{
00399         
00400     // Dimensions for the output arrays:
00401     int numFields = wedgeBasis.getCardinality();
00402     int numPoints = wedgeNodes.dimension(0);
00403     int spaceDim  = wedgeBasis.getBaseCellTopology().getDimension();
00404     
00405     // Generic array for values, grads, curls, etc. that will be properly sized before each call
00406     FieldContainer<double> vals;
00407     
00408     // Check VALUE of basis functions: resize vals to rank-2 container:
00409     vals.resize(numFields, numPoints);
00410     wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_VALUE);
00411     for (int i = 0; i < numFields; i++) {
00412       for (int j = 0; j < numPoints; j++) {
00413           int l =  i + j * numFields;
00414            if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
00415              errorFlag++;
00416              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00417 
00418              // Output the multi-index of the value where the error is:
00419              *outStream << " At multi-index { ";
00420              *outStream << i << " ";*outStream << j << " ";
00421              *outStream << "}  computed value: " << vals(i,j)
00422                << " but reference value: " << basisValues[l] << "\n";
00423          }
00424       }
00425     }
00426 
00427     
00428     
00429     // Check GRAD of basis function: resize vals to rank-3 container
00430     vals.resize(numFields, numPoints, spaceDim);
00431     wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_GRAD);
00432     for (int i = 0; i < numFields; i++) {
00433       for (int j = 0; j < numPoints; j++) {
00434         for (int k = 0; k < spaceDim; k++) {
00435           
00436           // basisGrads is (F,P,D), compute offset:
00437           int l = k + j * spaceDim + i * spaceDim * numPoints;
00438           if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00439             errorFlag++;
00440             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00441             
00442             // Output the multi-index of the value where the error is:
00443             *outStream << " At multi-index { ";
00444             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00445             *outStream << "}  computed grad component: " << vals(i,j,k)
00446               << " but reference grad component: " << basisGrads[l] << "\n";
00447           }
00448         }
00449       }
00450     }
00451     
00452     // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
00453     wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_D1);
00454     for (int i = 0; i < numFields; i++) {
00455       for (int j = 0; j < numPoints; j++) {
00456         for (int k = 0; k < spaceDim; k++) {
00457           
00458           // basisGrads is (F,P,D), compute offset:
00459           int l = k + j * spaceDim + i * spaceDim * numPoints;
00460           if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00461             errorFlag++;
00462             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00463             
00464             // Output the multi-index of the value where the error is:
00465             *outStream << " At multi-index { ";
00466             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00467             *outStream << "}  computed D1 component: " << vals(i,j,k)
00468               << " but reference D1 component: " << basisGrads[l] << "\n";
00469           }
00470         }
00471       }
00472     }
00473     
00474     
00475     // Check D2 of basis function
00476     int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim);
00477     vals.resize(numFields, numPoints, D2cardinality);    
00478     wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_D2);
00479     for (int i = 0; i < numFields; i++) {
00480       for (int j = 0; j < numPoints; j++) {
00481         for (int k = 0; k < D2cardinality; k++) {
00482           
00483           // basisD2 is (F,P,Dk), compute offset:
00484           int l = k + j * D2cardinality + i * D2cardinality * numPoints;
00485           if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) {
00486             errorFlag++;
00487             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00488             
00489             // Output the multi-index of the value where the error is:
00490             *outStream << " At multi-index { ";
00491             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00492             *outStream << "}  computed D2 component: " << vals(i,j,k)
00493               << " but reference D2 component: " << basisD2[l] << "\n";
00494           }
00495         }
00496       }
00497     }
00498     
00499     
00500     // Check D3 of basis function
00501     int D3cardinality = Intrepid::getDkCardinality(OPERATOR_D3, spaceDim);
00502     vals.resize(numFields, numPoints, D3cardinality);    
00503     wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_D3);
00504     
00505     for (int i = 0; i < numFields; i++) {
00506       for (int j = 0; j < numPoints; j++) {
00507         for (int k = 0; k < D3cardinality; k++) {
00508           
00509           // basisD3 is (F,P,Dk), compute offset:
00510           int l = k + j * D3cardinality + i * D3cardinality * numPoints;
00511           if (std::abs(vals(i,j,k) - basisD3[l]) > INTREPID_TOL) {
00512             errorFlag++;
00513             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00514             
00515             // Output the multi-index of the value where the error is:
00516             *outStream << " At multi-index { ";
00517             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00518             *outStream << "}  computed D3 component: " << vals(i,j,k)
00519               << " but reference D3 component: " << basisD3[l] << "\n";
00520           }
00521         }
00522       }
00523     }
00524     
00525     
00526     // Check D4 of basis function
00527     int D4cardinality = Intrepid::getDkCardinality(OPERATOR_D4, spaceDim);
00528     vals.resize(numFields, numPoints, D4cardinality);    
00529     wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_D4);
00530     for (int i = 0; i < numFields; i++) {
00531       for (int j = 0; j < numPoints; j++) {
00532         for (int k = 0; k < D4cardinality; k++) {
00533           
00534           // basisD4 is (F,P,Dk), compute offset:
00535           int l = k + j * D4cardinality + i * D4cardinality * numPoints;
00536           if (std::abs(vals(i,j,k) - basisD4[l]) > INTREPID_TOL) {
00537             errorFlag++;
00538             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00539             
00540             // Output the multi-index of the value where the error is:
00541             *outStream << " At multi-index { ";
00542             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00543             *outStream << "}  computed D4 component: " << vals(i,j,k)
00544               << " but reference D4 component: " << basisD2[l] << "\n";
00545           }
00546         }
00547       }
00548     }
00549     
00550         
00551     // Check all higher derivatives - must be zero. 
00552     for(EOperator op = OPERATOR_D5; op < OPERATOR_MAX; op++) {
00553       
00554       // The last dimension is the number of kth derivatives and needs to be resized for every Dk
00555       int DkCardin  = Intrepid::getDkCardinality(op, spaceDim);
00556       vals.resize(numFields, numPoints, DkCardin);    
00557 
00558       wedgeBasis.getValues(vals, wedgeNodes, op);
00559       for (int i = 0; i < vals.size(); i++) {
00560         if (std::abs(vals[i]) > INTREPID_TOL) {
00561           errorFlag++;
00562           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00563           
00564           // Get the multi-index of the value where the error is and the operator order
00565           std::vector<int> myIndex;
00566           vals.getMultiIndex(myIndex,i);
00567           int ord = Intrepid::getOperatorOrder(op);
00568           *outStream << " At multi-index { ";
00569           for(int j = 0; j < vals.rank(); j++) {
00570             *outStream << myIndex[j] << " ";
00571           }
00572           *outStream << "}  computed D"<< ord <<" component: " << vals[i] 
00573             << " but reference D" << ord << " component:  0 \n";
00574         }
00575       }
00576     }    
00577   }
00578   
00579   // Catch unexpected errors
00580   catch (std::logic_error err) {
00581     *outStream << err.what() << "\n\n";
00582     errorFlag = -1000;
00583   };
00584   
00585   if (errorFlag != 0)
00586     std::cout << "End Result: TEST FAILED\n";
00587   else
00588     std::cout << "End Result: TEST PASSED\n";
00589   
00590   // reset format state of std::cout
00591   std::cout.copyfmt(oldFormatState);
00592   
00593   return errorFlag;
00594 }