Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Basis/HGRAD_TET_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 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Pavel Bochev  (pbboche@sandia.gov),
00025 //                    Denis Ridzal  (dridzal@sandia.gov),
00026 //                    Kara Peterson (kjpeter@sandia.gov).
00027 //
00028 // ************************************************************************
00029 // @HEADER
00030 
00035 #include "Intrepid_FieldContainer.hpp"
00036 #include "Intrepid_HGRAD_TET_C2_FEM.hpp"
00037 #include "Teuchos_oblackholestream.hpp"
00038 #include "Teuchos_RCP.hpp"
00039 #include "Teuchos_GlobalMPISession.hpp"
00040 
00041 using namespace std;
00042 using namespace Intrepid;
00043 
00044 #define INTREPID_TEST_COMMAND( S , throwCounter, nException )                                                              \
00045 {                                                                                                                          \
00046   ++nException;                                                                                                            \
00047   try {                                                                                                                    \
00048     S ;                                                                                                                    \
00049   }                                                                                                                        \
00050   catch (std::logic_error err) {                                                                                           \
00051       ++throwCounter;                                                                                                      \
00052       *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
00053       *outStream << err.what() << '\n';                                                                                    \
00054       *outStream << "-------------------------------------------------------------------------------" << "\n\n";           \
00055   };                                                                                                                       \
00056 }
00057 
00058 int main(int argc, char *argv[]) {
00059   
00060   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00061 
00062   // This little trick lets us print to std::cout only if
00063   // a (dummy) command-line argument is provided.
00064   int iprint     = argc - 1;
00065   Teuchos::RCP<std::ostream> outStream;
00066   Teuchos::oblackholestream bhs; // outputs nothing
00067   if (iprint > 0)
00068     outStream = Teuchos::rcp(&std::cout, false);
00069   else
00070     outStream = Teuchos::rcp(&bhs, false);
00071   
00072   // Save the format state of the original std::cout.
00073   Teuchos::oblackholestream oldFormatState;
00074   oldFormatState.copyfmt(std::cout);
00075   
00076   *outStream \
00077     << "===============================================================================\n" \
00078     << "|                                                                             |\n" \
00079     << "|                 Unit Test (Basis_HGRAD_TET_C2_FEM)                          |\n" \
00080     << "|                                                                             |\n" \
00081     << "|     1) Conversion of Dof tags into Dof ordinals and back                    |\n" \
00082     << "|     2) Basis values for VALUE, GRAD, and Dk operators                       |\n" \
00083     << "|                                                                             |\n" \
00084     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00085     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00086     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00087     << "|                                                                             |\n" \
00088     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00089     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00090     << "|                                                                             |\n" \
00091     << "===============================================================================\n"\
00092     << "| TEST 1: Basis creation, exception testing                                   |\n"\
00093     << "===============================================================================\n";
00094   
00095   // Define basis and error flag
00096   Basis_HGRAD_TET_C2_FEM<double, FieldContainer<double> > tetBasis;
00097   int errorFlag = 0;
00098 
00099   // Initialize throw counter for exception testing
00100   int nException     = 0;
00101   int throwCounter   = 0;
00102 
00103   // Define array containing the 10 nodes of Tetrahedron<10>: 4 vertices + 6 edge midpoints  
00104   FieldContainer<double> tetNodes(10, 3);
00105   tetNodes(0,0) =  0.0;  tetNodes(0,1) =  0.0;  tetNodes(0,2) =  0.0;  
00106   tetNodes(1,0) =  1.0;  tetNodes(1,1) =  0.0;  tetNodes(1,2) =  0.0;  
00107   tetNodes(2,0) =  0.0;  tetNodes(2,1) =  1.0;  tetNodes(2,2) =  0.0;
00108   tetNodes(3,0) =  0.0;  tetNodes(3,1) =  0.0;  tetNodes(3,2) =  1.0;  
00109   
00110   tetNodes(4,0) =  0.5;  tetNodes(4,1) =  0.0;  tetNodes(4,2) =  0.0;
00111   tetNodes(5,0) =  0.5;  tetNodes(5,1) =  0.5;  tetNodes(5,2) =  0.0;  
00112   tetNodes(6,0) =  0.0;  tetNodes(6,1) =  0.5;  tetNodes(6,2) =  0.0;  
00113   tetNodes(7,0) =  0.0;  tetNodes(7,1) =  0.0;  tetNodes(7,2) =  0.5;  
00114   tetNodes(8,0) =  0.5;  tetNodes(8,1) =  0.0;  tetNodes(8,2) =  0.5;  
00115   tetNodes(9,0) =  0.0;  tetNodes(9,1) =  0.5;  tetNodes(9,2) =  0.5;  
00116 
00117 
00118   // Generic array for the output values; needs to be properly resized depending on the operator type
00119   FieldContainer<double> vals;
00120 
00121   try{
00122     // exception #1: CURL cannot be applied to scalar functions
00123     // resize vals to rank-3 container with dimensions (num. points, num. basis functions, arbitrary)
00124     vals.resize(tetBasis.getCardinality(), tetNodes.dimension(0), 4 );
00125     INTREPID_TEST_COMMAND( tetBasis.getValues(vals, tetNodes, OPERATOR_CURL), throwCounter, nException );
00126 
00127     // exception #2: DIV cannot be applied to scalar functions
00128     // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
00129     vals.resize(tetBasis.getCardinality(), tetNodes.dimension(0) );
00130     INTREPID_TEST_COMMAND( tetBasis.getValues(vals, tetNodes, OPERATOR_DIV), throwCounter, nException );
00131         
00132     // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 
00133     // getDofTag() to access invalid array elements thereby causing bounds check exception
00134     // exception #3
00135     INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(3,0,0), throwCounter, nException );
00136     // exception #4
00137     INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(1,1,1), throwCounter, nException );
00138     // exception #5
00139     INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(0,4,0), throwCounter, nException );
00140     // exception #6
00141     INTREPID_TEST_COMMAND( tetBasis.getDofTag(10), throwCounter, nException );
00142     // exception #7
00143     INTREPID_TEST_COMMAND( tetBasis.getDofTag(-1), throwCounter, nException );
00144     
00145 #ifdef HAVE_INTREPID_DEBUG 
00146     // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays
00147     // exception #8: input points array must be of rank-2
00148     FieldContainer<double> badPoints1(4, 5, 3);
00149     INTREPID_TEST_COMMAND( tetBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00150     
00151     // exception #9 dimension 1 in the input point array must equal space dimension of the cell
00152     FieldContainer<double> badPoints2(4, tetBasis.getBaseCellTopology().getDimension() - 1);
00153     INTREPID_TEST_COMMAND( tetBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00154     
00155     // exception #10 output values must be of rank-2 for OPERATOR_VALUE
00156     FieldContainer<double> badVals1(4, 3, 1);
00157     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals1, tetNodes, OPERATOR_VALUE), throwCounter, nException );
00158     
00159     // exception #11 output values must be of rank-3 for OPERATOR_GRAD
00160     FieldContainer<double> badVals2(4, 3);
00161     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals2, tetNodes, OPERATOR_GRAD), throwCounter, nException );
00162     
00163     // exception #12 output values must be of rank-3 for OPERATOR_D1
00164     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals2, tetNodes, OPERATOR_D1), throwCounter, nException );
00165     
00166     // exception #13 output values must be of rank-3 for OPERATOR_D2
00167     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals2, tetNodes, OPERATOR_D2), throwCounter, nException );
00168     
00169     // exception #14 incorrect 0th dimension of output array (must equal number of basis functions)
00170     FieldContainer<double> badVals3(tetBasis.getCardinality() + 1, tetNodes.dimension(0));
00171     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals3, tetNodes, OPERATOR_VALUE), throwCounter, nException );
00172     
00173     // exception #15 incorrect 1st dimension of output array (must equal number of points)
00174     FieldContainer<double> badVals4(tetBasis.getCardinality(), tetNodes.dimension(0) + 1);
00175     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals4, tetNodes, OPERATOR_VALUE), throwCounter, nException );
00176     
00177     // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
00178     FieldContainer<double> badVals5(tetBasis.getCardinality(), tetNodes.dimension(0), tetBasis.getBaseCellTopology().getDimension() + 1);
00179     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals5, tetNodes, OPERATOR_GRAD), throwCounter, nException );
00180     
00181     // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 2D)
00182     FieldContainer<double> badVals6(tetBasis.getCardinality(), tetNodes.dimension(0), 40);
00183     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals6, tetNodes, OPERATOR_D1), throwCounter, nException );
00184     
00185     // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 2D)
00186     INTREPID_TEST_COMMAND( tetBasis.getValues(badVals6, tetNodes, OPERATOR_D2), throwCounter, nException );
00187 #endif
00188     
00189   }
00190   catch (std::logic_error err) {
00191     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00192     *outStream << err.what() << '\n';
00193     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00194     errorFlag = -1000;
00195   };
00196   
00197   // Check if number of thrown exceptions matches the one we expect 
00198   if (throwCounter != nException) {
00199     errorFlag++;
00200     *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00201   }
00202   
00203   *outStream \
00204     << "\n"
00205     << "===============================================================================\n"\
00206     << "| TEST 2: correctness of tag to enum and enum to tag lookups                  |\n"\
00207     << "===============================================================================\n";
00208   
00209   try{
00210     std::vector<std::vector<int> > allTags = tetBasis.getAllDofTags();
00211     
00212     // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
00213     for (unsigned i = 0; i < allTags.size(); i++) {
00214       int bfOrd  = tetBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
00215       
00216       std::vector<int> myTag = tetBasis.getDofTag(bfOrd);
00217        if( !( (myTag[0] == allTags[i][0]) &&
00218               (myTag[1] == allTags[i][1]) &&
00219               (myTag[2] == allTags[i][2]) &&
00220               (myTag[3] == allTags[i][3]) ) ) {
00221         errorFlag++;
00222         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00223         *outStream << " getDofOrdinal( {" 
00224           << allTags[i][0] << ", " 
00225           << allTags[i][1] << ", " 
00226           << allTags[i][2] << ", " 
00227           << allTags[i][3] << "}) = " << bfOrd <<" but \n";   
00228         *outStream << " getDofTag(" << bfOrd << ") = { "
00229           << myTag[0] << ", " 
00230           << myTag[1] << ", " 
00231           << myTag[2] << ", " 
00232           << myTag[3] << "}\n";        
00233       }
00234     }
00235     
00236     // Now do the same but loop over basis functions
00237     for( int bfOrd = 0; bfOrd < tetBasis.getCardinality(); bfOrd++) {
00238       std::vector<int> myTag  = tetBasis.getDofTag(bfOrd);
00239       int myBfOrd = tetBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
00240       if( bfOrd != myBfOrd) {
00241         errorFlag++;
00242         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00243         *outStream << " getDofTag(" << bfOrd << ") = { "
00244           << myTag[0] << ", " 
00245           << myTag[1] << ", " 
00246           << myTag[2] << ", " 
00247           << myTag[3] << "} but getDofOrdinal({" 
00248           << myTag[0] << ", " 
00249           << myTag[1] << ", " 
00250           << myTag[2] << ", " 
00251           << myTag[3] << "} ) = " << myBfOrd << "\n";
00252       }
00253     }
00254   }
00255   catch (std::logic_error err){
00256     *outStream << err.what() << "\n\n";
00257     errorFlag = -1000;
00258   };
00259   
00260   *outStream \
00261     << "\n"
00262     << "===============================================================================\n"\
00263     << "| TEST 3: correctness of basis function values                                |\n"\
00264     << "===============================================================================\n";
00265   
00266   outStream -> precision(20);
00267   
00268   // VALUE: in (F,P) format
00269   double basisValues[] = {
00270     1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, \
00271     0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, \
00272     0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, \
00273     0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00274     1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, \
00275     0, 0, 0, 1.00000 };
00276   
00277   // GRAD and D1: in (F,P,D) format
00278   double basisGrads[] = {
00279     -3.00000, -3.00000, -3.00000, 1.00000, 1.00000, 1.00000, 1.00000, \
00280     1.00000, 1.00000, 1.00000, 1.00000, 1.00000, -1.00000, -1.00000, \
00281     -1.00000, 1.00000, 1.00000, 1.00000, -1.00000, -1.00000, -1.00000, \
00282     -1.00000, -1.00000, -1.00000, 1.00000, 1.00000, 1.00000, 1.00000, \
00283     1.00000, 1.00000, -1.00000, 0, 0, 3.00000, 0, 0, -1.00000, 0, 0, \
00284     -1.00000, 0, 0, 1.00000, 0, 0, 1.00000, 0, 0, -1.00000, 0, 0, \
00285     -1.00000, 0, 0, 1.00000, 0, 0, -1.00000, 0, 0, 0, -1.00000, 0, 0, \
00286     -1.00000, 0, 0, 3.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \
00287     1.00000, 0, 0, 1.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \
00288     1.00000, 0, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \
00289     3.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \
00290     1.00000, 0, 0, 1.00000, 0, 0, 1.00000, 4.00000, 0, 0, -4.00000, \
00291     -4.00000, -4.00000, 0, 0, 0, 0, 0, 0, 0, -2.00000, -2.00000, \
00292     -2.00000, -2.00000, -2.00000, 2.00000, 0, 0, 2.00000, 0, 0, -2.00000, \
00293     -2.00000, -2.00000, 0, 0, 0, 0, 0, 0, 0, 4.00000, 0, 4.00000, 0, 0, \
00294     0, 0, 0, 0, 2.00000, 0, 2.00000, 2.00000, 0, 2.00000, 0, 0, 0, 0, 0, \
00295     0, 2.00000, 0, 2.00000, 0, 0, 0, 4.00000, 0, 0, 0, 0, -4.00000, \
00296     -4.00000, -4.00000, 0, 0, 0, 0, 2.00000, 0, -2.00000, -2.00000, \
00297     -2.00000, -2.00000, 0, -2.00000, 0, 2.00000, 0, 0, 0, 0, -2.00000, \
00298     -2.00000, -2.00000, 0, 0, 4.00000, 0, 0, 0, 0, 0, 0, -4.00000, \
00299     -4.00000, -4.00000, 0, 0, 2.00000, 0, 0, 0, 0, 0, 2.00000, -2.00000, \
00300     -2.00000, 0, -2.00000, -2.00000, -2.00000, -2.00000, -2.00000, \
00301     -2.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 4.00000, 0, 0, 0, 0, \
00302     2.00000, 0, 0, 2.00000, 0, 0, 0, 2.00000, 0, 0, 2.00000, 0, 2.00000, \
00303     2.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 4.00000, 0, 4.00000, 0, 0, 0, \
00304     0, 0, 0, 2.00000, 0, 0, 2.00000, 0, 2.00000, 0, 0, 2.00000, 0, 0, \
00305     2.00000, 2.00000};
00306   
00307   // D2 values in (F,P, Dk) format
00308   double basisD2[]={
00309     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00310     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00311     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00312     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00313     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00314     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00315     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00316     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \
00317     4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 0, 0, 0, 0, 0, 4.00000, \
00318     0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, \
00319     4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \
00320     0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \
00321     0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, \
00322     4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \
00323     0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \
00324     0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 0, 0, 4.00000, \
00325     0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, \
00326     4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \
00327     0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \
00328     0, 0, 4.00000, -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, \
00329     -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, \
00330     -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, \
00331     0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, \
00332     -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, \
00333     -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, \
00334     0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \
00335     0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, \
00336     0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, \
00337     0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, -4.00000, 0, -8.00000, -4.00000, \
00338     0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, -8.00000, \
00339     -4.00000, 0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, \
00340     -8.00000, -4.00000, 0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, \
00341     -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, -8.00000, \
00342     -4.00000, 0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, \
00343     -8.00000, -4.00000, 0, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, \
00344     -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, -4.00000, \
00345     -8.00000, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, \
00346     -4.00000, -8.00000, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, \
00347     -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, -4.00000, \
00348     -8.00000, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, \
00349     -4.00000, -8.00000, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \
00350     0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, \
00351     0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, \
00352     0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 0, 0, \
00353     4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \
00354     0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \
00355     0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, \
00356     0, 0, 0, 4.00000, 0
00357   };
00358   
00359   
00360   try{
00361         
00362     // Dimensions for the output arrays:
00363     int numFields = tetBasis.getCardinality();
00364     int numPoints = tetNodes.dimension(0);
00365     int spaceDim  = tetBasis.getBaseCellTopology().getDimension();
00366     
00367     // Generic array for values, grads, curls, etc. that will be properly sized before each call
00368     FieldContainer<double> vals;
00369     
00370     // Check VALUE of basis functions: resize vals to rank-2 container:
00371     vals.resize(numFields, numPoints);
00372     tetBasis.getValues(vals, tetNodes, OPERATOR_VALUE);
00373     for (int i = 0; i < numFields; i++) {
00374       for (int j = 0; j < numPoints; j++) {
00375           int l =  i + j * numFields;
00376            if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
00377              errorFlag++;
00378              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00379 
00380              // Output the multi-index of the value where the error is:
00381              *outStream << " At multi-index { ";
00382              *outStream << i << " ";*outStream << j << " ";
00383              *outStream << "}  computed value: " << vals(i,j)
00384                << " but reference value: " << basisValues[l] << "\n";
00385          }
00386       }
00387     }
00388 
00389     // Check GRAD of basis function: resize vals to rank-3 container
00390     vals.resize(numFields, numPoints, spaceDim);
00391     tetBasis.getValues(vals, tetNodes, OPERATOR_GRAD);
00392     for (int i = 0; i < numFields; i++) {
00393       for (int j = 0; j < numPoints; j++) {
00394         for (int k = 0; k < spaceDim; k++) {
00395  
00396           // basisGrads is (F,P,D), compute offset:
00397           int l = k + j * spaceDim + i * spaceDim * numPoints;
00398            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00399              errorFlag++;
00400              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00401 
00402              // Output the multi-index of the value where the error is:
00403              *outStream << " At multi-index { ";
00404              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00405              *outStream << "}  computed grad component: " << vals(i,j,k)
00406                << " but reference grad component: " << basisGrads[l] << "\n";
00407             }
00408          }
00409       }
00410     }
00411 
00412     // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
00413     tetBasis.getValues(vals, tetNodes, OPERATOR_D1);
00414     for (int i = 0; i < numFields; i++) {
00415       for (int j = 0; j < numPoints; j++) {
00416         for (int k = 0; k < spaceDim; k++) {
00417           
00418           // basisGrads is (F,P,D), compute offset:
00419           int l = k + j * spaceDim + i * spaceDim * numPoints;
00420            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00421              errorFlag++;
00422              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00423 
00424              // Output the multi-index of the value where the error is:
00425              *outStream << " At multi-index { ";
00426              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00427              *outStream << "}  computed D1 component: " << vals(i,j,k)
00428                << " but reference D1 component: " << basisGrads[l] << "\n";
00429             }
00430          }
00431       }
00432     }
00433 
00434     // Check D2 of basis function
00435     int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim);
00436     vals.resize(numFields, numPoints, D2cardinality);    
00437     tetBasis.getValues(vals, tetNodes, OPERATOR_D2);
00438     for (int i = 0; i < numFields; i++) {
00439       for (int j = 0; j < numPoints; j++) {
00440         for (int k = 0; k < D2cardinality; k++) {
00441           
00442           // basisD2 is (F,P,Dk), compute offset:
00443           int l = k + j * D2cardinality + i * D2cardinality * numPoints;
00444           if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) {
00445             errorFlag++;
00446             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00447             
00448             // Output the multi-index of the value where the error is:
00449             *outStream << " At multi-index { ";
00450             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00451             *outStream << "}  computed D2 component: " << vals(i,j,k)
00452               << " but reference D2 component: " << basisD2[l] << "\n";
00453           }
00454         }
00455       }
00456     }
00457     
00458     
00459     // Check all higher derivatives - must be zero. 
00460     for(EOperator op = OPERATOR_D3; op < OPERATOR_MAX; op++) {
00461       
00462       // The last dimension is the number of kth derivatives and needs to be resized for every Dk
00463       int DkCardin  = Intrepid::getDkCardinality(op, spaceDim);
00464       vals.resize(numFields, numPoints, DkCardin);    
00465 
00466       tetBasis.getValues(vals, tetNodes, op);
00467       for (int i = 0; i < vals.size(); i++) {
00468         if (std::abs(vals[i]) > INTREPID_TOL) {
00469           errorFlag++;
00470           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00471           
00472           // Get the multi-index of the value where the error is and the operator order
00473           std::vector<int> myIndex;
00474           vals.getMultiIndex(myIndex,i);
00475           int ord = Intrepid::getOperatorOrder(op);
00476           *outStream << " At multi-index { ";
00477           for(int j = 0; j < vals.rank(); j++) {
00478             *outStream << myIndex[j] << " ";
00479           }
00480           *outStream << "}  computed D"<< ord <<" component: " << vals[i] 
00481             << " but reference D" << ord << " component:  0 \n";
00482         }
00483       }
00484     }    
00485   }
00486   
00487   // Catch unexpected errors
00488   catch (std::logic_error err) {
00489     *outStream << err.what() << "\n\n";
00490     errorFlag = -1000;
00491   };
00492   
00493   if (errorFlag != 0)
00494     std::cout << "End Result: TEST FAILED\n";
00495   else
00496     std::cout << "End Result: TEST PASSED\n";
00497   
00498   // reset format state of std::cout
00499   std::cout.copyfmt(oldFormatState);
00500   
00501   return errorFlag;
00502 }