Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Basis/HGRAD_HEX_Cn_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_HEX_Cn_FEM.hpp"
00037 #include "Intrepid_PointTools.hpp"
00038 #include "Teuchos_oblackholestream.hpp"
00039 #include "Teuchos_RCP.hpp"
00040 #include "Teuchos_GlobalMPISession.hpp"
00041 
00042 using namespace std;
00043 using namespace Intrepid;
00044 
00045 #define INTREPID_TEST_COMMAND( S , throwCounter, nException )                                                              \
00046 {                                                                                                                          \
00047   ++nException;                                                                                                            \
00048   try {                                                                                                                    \
00049     S ;                                                                                                                    \
00050   }                                                                                                                        \
00051   catch (std::logic_error err) {                                                                                           \
00052       ++throwCounter;                                                                                                      \
00053       *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
00054       *outStream << err.what() << '\n';                                                                                    \
00055       *outStream << "-------------------------------------------------------------------------------" << "\n\n";           \
00056   };                                                                                                                       \
00057 }
00058 
00059 int main(int argc, char *argv[]) {
00060   
00061   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00062 
00063   // This little trick lets us print to std::cout only if
00064   // a (dummy) command-line argument is provided.
00065   int iprint     = argc - 1;
00066   Teuchos::RCP<std::ostream> outStream;
00067   Teuchos::oblackholestream bhs; // outputs nothing
00068   if (iprint > 0)
00069     outStream = Teuchos::rcp(&std::cout, false);
00070   else
00071     outStream = Teuchos::rcp(&bhs, false);
00072   
00073   // Save the format state of the original std::cout.
00074   Teuchos::oblackholestream oldFormatState;
00075   oldFormatState.copyfmt(std::cout);
00076   
00077   *outStream \
00078     << "===============================================================================\n" \
00079     << "|                                                                             |\n" \
00080     << "|                 Unit Test (Basis_HGRAD_HEX_C2_FEM)                          |\n" \
00081     << "|                                                                             |\n" \
00082     << "|     1) Conversion of Dof tags into Dof ordinals and back                    |\n" \
00083     << "|     2) Basis values for VALUE, GRAD, and Dk operators                       |\n" \
00084     << "|                                                                             |\n" \
00085     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00086     << "|                      Robert Kirby  (robert.c.kirby@ttu.edu),                |\n" \
00087     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00088     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00089     << "|                                                                             |\n" \
00090     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00091     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00092     << "|                                                                             |\n" \
00093     << "===============================================================================\n"\
00094     << "| TEST 1: Basis creation, exception testing                                   |\n"\
00095     << "===============================================================================\n";
00096   
00097   // Define basis and error flag
00098   // get points
00099   const int deg=2;
00100   shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >()); 
00101   FieldContainer<double> pts(PointTools::getLatticeSize(line,deg),1);
00102   PointTools::getLattice<double,FieldContainer<double> >(pts,line,deg);
00103 
00104   Basis_HGRAD_HEX_Cn_FEM<double, FieldContainer<double> > hexBasis(deg,deg,deg,pts,pts,pts);
00105 
00106   int errorFlag = 0;
00107 
00108   // Initialize throw counter for exception testing
00109   int nException     = 0;
00110   int throwCounter   = 0;
00111 
00112   // Define arrayS containing the 27 nodes of hexahedron<27> topology
00113   FieldContainer<double> hexNodes(27, 3);
00114   // do it lexicographically as a lattice
00115   hexNodes(0, 0) = -1.0;   hexNodes(0, 1) = -1.0;  hexNodes(0, 2) = -1.0;
00116   hexNodes(1, 0) =  0.0;   hexNodes(1, 1) = -1.0;  hexNodes(1, 2) = -1.0;
00117   hexNodes(2, 0) =  1.0;   hexNodes(2, 1) = -1.0;  hexNodes(2, 2) = -1.0;
00118   hexNodes(3, 0) = -1.0;   hexNodes(3, 1) =  0.0;  hexNodes(3, 2) = -1.0;
00119   hexNodes(4, 0) =  0.0;   hexNodes(4, 1) =  0.0;  hexNodes(4, 2) = -1.0;
00120   hexNodes(5, 0) =  1.0;   hexNodes(5, 1) =  0.0;  hexNodes(5, 2) = -1.0;
00121   hexNodes(6, 0) = -1.0;   hexNodes(6, 1) =  1.0;  hexNodes(6, 2) = -1.0;
00122   hexNodes(7, 0) = 0.0;    hexNodes(7, 1) =  1.0;  hexNodes(7, 2) = -1.0;
00123   hexNodes(8, 0) = 1.0;    hexNodes(8, 1) =  1.0;  hexNodes(8, 2) = -1.0;
00124   hexNodes(9, 0) = -1.0;   hexNodes(9, 1) = -1.0;  hexNodes(9, 2) = 0.0;
00125   hexNodes(10, 0) =  0.0;   hexNodes(10, 1) = -1.0;  hexNodes(10, 2) = 0.0;
00126   hexNodes(11, 0) =  1.0;   hexNodes(11, 1) = -1.0;  hexNodes(11, 2) = 0.0;
00127   hexNodes(12, 0) = -1.0;   hexNodes(12, 1) =  0.0;  hexNodes(12, 2) = 0.0;
00128   hexNodes(13, 0) =  0.0;   hexNodes(13, 1) =  0.0;  hexNodes(13, 2) = 0.0;
00129   hexNodes(14, 0) =  1.0;   hexNodes(14, 1) =  0.0;  hexNodes(14, 2) = 0.0;
00130   hexNodes(15, 0) = -1.0;   hexNodes(15, 1) =  1.0;  hexNodes(15, 2) = 0.0;
00131   hexNodes(16, 0) = 0.0;    hexNodes(16, 1) =  1.0;  hexNodes(16, 2) = 0.0;
00132   hexNodes(17, 0) = 1.0;    hexNodes(17, 1) =  1.0;  hexNodes(17, 2) = 0.0;
00133   hexNodes(18, 0) = -1.0;   hexNodes(18, 1) = -1.0;  hexNodes(18, 2) = 1.0;
00134   hexNodes(19, 0) =  0.0;   hexNodes(19, 1) = -1.0;  hexNodes(19, 2) = 1.0;
00135   hexNodes(20, 0) =  1.0;   hexNodes(20, 1) = -1.0;  hexNodes(20, 2) = 1.0;
00136   hexNodes(21, 0) = -1.0;   hexNodes(21, 1) =  0.0;  hexNodes(21, 2) = 1.0;
00137   hexNodes(22, 0) =  0.0;   hexNodes(22, 1) =  0.0;  hexNodes(22, 2) = 1.0;
00138   hexNodes(23, 0) =  1.0;   hexNodes(23, 1) =  0.0;  hexNodes(23, 2) = 1.0;
00139   hexNodes(24, 0) = -1.0;   hexNodes(24, 1) =  1.0;  hexNodes(24, 2) = 1.0;
00140   hexNodes(25, 0) = 0.0;    hexNodes(25, 1) =  1.0;  hexNodes(25, 2) = 1.0;
00141   hexNodes(26, 0) = 1.0;    hexNodes(26, 1) =  1.0;  hexNodes(26, 2) = 1.0;
00142 
00143  
00144   // Generic array for the output values; needs to be properly resized depending on the operator type
00145   FieldContainer<double> vals;
00146 
00147   try{
00148     // exception #1: CURL cannot be applied to scalar functions in 3D
00149     // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
00150     vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0), 4 );
00151     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_CURL), throwCounter, nException );
00152 
00153     // exception #2: DIV cannot be applied to scalar functions in 3D
00154     // resize vals to rank-2 container with dimensions (num. basis functions, num. points)
00155     vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0) );
00156     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_DIV), throwCounter, nException );
00157         
00158     // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 
00159     // getDofTag() to access invalid array elements thereby causing bounds check exception
00160     // exception #3
00161     INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(3,10,0), throwCounter, nException );    
00162     // exception #4
00163     INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(1,2,1), throwCounter, nException );
00164     // exception #5
00165     INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(0,4,1), throwCounter, nException );
00166     // exception #6
00167     INTREPID_TEST_COMMAND( hexBasis.getDofTag(28), throwCounter, nException );
00168     // exception #7
00169     INTREPID_TEST_COMMAND( hexBasis.getDofTag(-1), throwCounter, nException );
00170 
00171 #ifdef HAVE_INTREPID_DEBUG 
00172     // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays
00173     // exception #8: input points array must be of rank-2
00174     FieldContainer<double> badPoints1(4, 5, 3);
00175     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00176 
00177     // exception #9 dimension 1 in the input point array must equal space dimension of the cell
00178     FieldContainer<double> badPoints2(4, hexBasis.getBaseCellTopology().getDimension() - 1);
00179     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00180 
00181     // exception #10 output values must be of rank-2 for OPERATOR_VALUE
00182     FieldContainer<double> badVals1(4, 3, 1);
00183     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_VALUE), throwCounter, nException );
00184  
00185     // exception #11 output values must be of rank-3 for OPERATOR_GRAD
00186     FieldContainer<double> badVals2(4, 3);
00187     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_GRAD), throwCounter, nException );
00188     
00189     // exception #12 output values must be of rank-3 for OPERATOR_D1
00190     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D1), throwCounter, nException );
00191     
00192     // exception #13 output values must be of rank-3 for OPERATOR_D2
00193     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D2), throwCounter, nException );
00194     
00195     // exception #14 incorrect 0th dimension of output array (must equal number of basis functions)
00196     FieldContainer<double> badVals3(hexBasis.getCardinality() + 1, hexNodes.dimension(0));
00197     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals3, hexNodes, OPERATOR_VALUE), throwCounter, nException );
00198     
00199     // exception #15 incorrect 1st dimension of output array (must equal number of points)
00200     FieldContainer<double> badVals4(hexBasis.getCardinality(), hexNodes.dimension(0) + 1);
00201     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_VALUE), throwCounter, nException );
00202     
00203     // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
00204     FieldContainer<double> badVals5(hexBasis.getCardinality(), hexNodes.dimension(0), hexBasis.getBaseCellTopology().getDimension() - 1);
00205     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals5, hexNodes, OPERATOR_GRAD), throwCounter, nException );
00206     
00207     // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 3D)
00208     FieldContainer<double> badVals6(hexBasis.getCardinality(), hexNodes.dimension(0), 40);
00209     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals6, hexNodes, OPERATOR_D2), throwCounter, nException );
00210     
00211     // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 3D)
00212     FieldContainer<double> badVals7(hexBasis.getCardinality(), hexNodes.dimension(0), 50);
00213     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals7, hexNodes, OPERATOR_D3), throwCounter, nException );
00214 #endif
00215 
00216   }
00217   catch (std::logic_error err) {
00218     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00219     *outStream << err.what() << '\n';
00220     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00221     errorFlag = -1000;
00222   };
00223   
00224   // Check if number of thrown exceptions matches the one we expect 
00225   // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match.
00226   if (throwCounter != nException) {
00227     errorFlag++;
00228     *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00229   }
00230   
00231   *outStream \
00232     << "\n"
00233     << "===============================================================================\n"\
00234     << "| TEST 2: correctness of tag to enum and enum to tag lookups                  |\n"\
00235     << "===============================================================================\n";
00236   
00237   try{
00238     std::vector<std::vector<int> > allTags = hexBasis.getAllDofTags();
00239     
00240     // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
00241     for (unsigned i = 0; i < allTags.size(); i++) {
00242       int bfOrd  = hexBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
00243       
00244 
00245 
00246       std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
00247        if( !( (myTag[0] == allTags[i][0]) &&
00248               (myTag[1] == allTags[i][1]) &&
00249               (myTag[2] == allTags[i][2]) &&
00250               (myTag[3] == allTags[i][3]) ) ) {
00251         errorFlag++;
00252         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00253         *outStream << " getDofOrdinal( {" 
00254           << allTags[i][0] << ", " 
00255           << allTags[i][1] << ", " 
00256           << allTags[i][2] << ", " 
00257           << allTags[i][3] << "}) = " << bfOrd <<" but \n";   
00258         *outStream << " getDofTag(" << bfOrd << ") = { "
00259           << myTag[0] << ", " 
00260           << myTag[1] << ", " 
00261           << myTag[2] << ", " 
00262           << myTag[3] << "}\n";        
00263       }
00264     }
00265     
00266     // Now do the same but loop over basis functions
00267     for( int bfOrd = 0; bfOrd < hexBasis.getCardinality(); bfOrd++) {
00268       std::vector<int> myTag  = hexBasis.getDofTag(bfOrd);
00269       int myBfOrd = hexBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
00270       if( bfOrd != myBfOrd) {
00271         errorFlag++;
00272         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00273         *outStream << " getDofTag(" << bfOrd << ") = { "
00274           << myTag[0] << ", " 
00275           << myTag[1] << ", " 
00276           << myTag[2] << ", " 
00277           << myTag[3] << "} but getDofOrdinal({" 
00278           << myTag[0] << ", " 
00279           << myTag[1] << ", " 
00280           << myTag[2] << ", " 
00281           << myTag[3] << "} ) = " << myBfOrd << "\n";
00282       }
00283     }
00284   }
00285   catch (std::logic_error err){
00286     *outStream << err.what() << "\n\n";
00287     errorFlag = -1000;
00288   };
00289   
00290   
00291   *outStream \
00292     << "\n"
00293     << "===============================================================================\n"\
00294     << "| TEST 3: correctness of basis function values                                |\n"\
00295     << "===============================================================================\n";
00296   
00297   outStream -> precision(20);
00298   
00299   // VALUE: Each row gives the 27 correct basis set values at an evaluation point
00300   double basisValues[] = {
00301     1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00302     0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00303     0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00304     0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00305     0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00306     0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00307     0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00308     0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00309     0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00310     0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00311     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00312     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00313     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00314     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00315     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00316     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00317     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00318     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00319     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, \
00320     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, \
00321     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, \
00322     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, \
00323     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, \
00324     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, \
00325     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, \
00326     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, \
00327     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000  };
00328   
00329   
00330   // GRAD, D1, D2, D3 and D4 test values are stored in files due to their large size
00331   std::string     fileName;
00332   std::ifstream   dataFile;
00333 
00334   // GRAD and D1 values are stored in (F,P,D) format in a data file. Read file and do the test
00335   std::vector<double> basisGrads;           // Flat array for the gradient values.
00336   
00337   fileName = "./testdata/HEX_C2_GradVals.dat";
00338   dataFile.open(fileName.c_str());
00339   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00340                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open GRAD values data file, test aborted.");
00341   while (!dataFile.eof() ){
00342     double temp;
00343     string line;                            // string for one line of input file
00344     std::getline(dataFile, line);           // get next line from file
00345     stringstream data_line(line);           // convert to stringstream
00346     while(data_line >> temp){               // extract value from line
00347       basisGrads.push_back(temp);           // push into vector
00348     }
00349   }
00350   // It turns out that just closing and then opening the ifstream variable does not reset it
00351   // and subsequent open() command fails. One fix is to explicitely clear the ifstream, or
00352   // scope the variables.
00353   dataFile.close();
00354   dataFile.clear();
00355     
00356    
00357   //D2: flat array with the values of D2 applied to basis functions. Multi-index is (F,P,D2cardinality)
00358   std::vector<double> basisD2; 
00359   fileName = "./testdata/HEX_C2_D2Vals.dat";  
00360   dataFile.open(fileName.c_str());
00361   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00362                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open D2 values data file, test aborted.");
00363   while (!dataFile.eof() ){
00364     double temp;
00365     string line;                            // string for one line of input file
00366     std::getline(dataFile, line);           // get next line from file
00367     stringstream data_line(line);           // convert to stringstream
00368     while(data_line >> temp){               // extract value from line
00369       basisD2.push_back(temp);              // push into vector
00370     }
00371   }
00372   dataFile.close();
00373   dataFile.clear();
00374 
00375   
00376   //D3: flat array with the values of D3 applied to basis functions. Multi-index is (F,P,D3cardinality)
00377   std::vector<double> basisD3;
00378   
00379   fileName = "./testdata/HEX_C2_D3Vals.dat";  
00380   dataFile.open(fileName.c_str());
00381   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00382                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open D3 values data file, test aborted.");
00383   
00384   while (!dataFile.eof() ){
00385     double temp;
00386     string line;                            // string for one line of input file
00387     std::getline(dataFile, line);           // get next line from file
00388     stringstream data_line(line);           // convert to stringstream
00389     while(data_line >> temp){               // extract value from line
00390       basisD3.push_back(temp);              // push into vector
00391     }
00392   }
00393   dataFile.close();
00394   dataFile.clear();
00395  
00396   
00397   //D4: flat array with the values of D applied to basis functions. Multi-index is (F,P,D4cardinality)
00398   std::vector<double> basisD4;
00399   
00400   fileName = "./testdata/HEX_C2_D4Vals.dat";
00401   dataFile.open(fileName.c_str());
00402   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00403                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open D4 values data file, test aborted.");
00404   
00405   while (!dataFile.eof() ){
00406     double temp;
00407     string line;                            // string for one line of input file
00408     std::getline(dataFile, line);           // get next line from file
00409     stringstream data_line(line);           // convert to stringstream
00410     while(data_line >> temp){               // extract value from line
00411       basisD4.push_back(temp);              // push into vector
00412     }
00413   }
00414   dataFile.close();
00415   dataFile.clear();
00416   
00417 
00418   try{
00419         
00420     // Dimensions for the output arrays:
00421     int numFields = hexBasis.getCardinality();
00422     int numPoints = hexNodes.dimension(0);
00423     int spaceDim  = hexBasis.getBaseCellTopology().getDimension();
00424     
00425     // Generic array for values, grads, curls, etc. that will be properly sized before each call
00426     FieldContainer<double> vals;
00427     
00428     // Check VALUE of basis functions: resize vals to rank-2 container:
00429     vals.resize(numFields, numPoints);
00430     hexBasis.getValues(vals, hexNodes, OPERATOR_VALUE);
00431     for (int i = 0; i < numFields; i++) {
00432       for (int j = 0; j < numPoints; j++) {
00433         int l =  i + j * numFields;
00434         if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
00435           errorFlag++;
00436           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00437           
00438           // Output the multi-index of the value where the error is:
00439           *outStream << " At multi-index { ";
00440           *outStream << i << " ";*outStream << j << " ";
00441           *outStream << "}  computed value: " << vals(i,j)
00442                      << " but reference value: " << basisValues[l] << "\n";
00443         }
00444       }
00445     }
00446     
00447     // Check GRAD of basis function: resize vals to rank-3 container
00448     vals.resize(numFields, numPoints, spaceDim);
00449     hexBasis.getValues(vals, hexNodes, OPERATOR_GRAD);
00450     for (int i = 0; i < numFields; i++) {
00451       for (int j = 0; j < numPoints; j++) {
00452         for (int k = 0; k < spaceDim; k++) {
00453           
00454           // basisGrads is (F,P,D), compute offset:
00455           int l = k + j * spaceDim + i * spaceDim * numPoints;
00456            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00457              errorFlag++;
00458              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00459 
00460              // Output the multi-index of the value where the error is:
00461              *outStream << " At multi-index { ";
00462              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00463              *outStream << "}  computed grad component: " << vals(i,j,k)
00464                << " but reference grad component: " << basisGrads[l] << "\n";
00465             }
00466          }
00467       }
00468     }
00469     
00470     // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
00471     hexBasis.getValues(vals, hexNodes, OPERATOR_D1);
00472     for (int i = 0; i < numFields; i++) {
00473       for (int j = 0; j < numPoints; j++) {
00474         for (int k = 0; k < spaceDim; k++) {
00475 
00476           // basisGrads is (F,P,D), compute offset:
00477           int l = k + j * spaceDim + i * spaceDim * numPoints;
00478            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00479              errorFlag++;
00480              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00481 
00482              // Output the multi-index of the value where the error is:
00483              *outStream << " At multi-index { ";
00484              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00485              *outStream << "}  computed D1 component: " << vals(i,j,k)
00486                << " but reference D1 component: " << basisGrads[l] << "\n";
00487             }
00488          }
00489       }
00490     }
00491 
00492     
00493     // Check D2 of basis function
00494     int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim);
00495     vals.resize(numFields, numPoints, D2cardinality);    
00496     hexBasis.getValues(vals, hexNodes, OPERATOR_D2);
00497     for (int i = 0; i < numFields; i++) {
00498       for (int j = 0; j < numPoints; j++) {
00499         for (int k = 0; k < D2cardinality; k++) {
00500  
00501           // basisD2 is (F,P,Dk), compute offset:
00502           int l = k + j * D2cardinality + i * D2cardinality * numPoints;
00503            if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) {
00504              errorFlag++;
00505              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00506 
00507              // Output the multi-index of the value where the error is:
00508              *outStream << " At multi-index { ";
00509              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00510              *outStream << "}  computed D2 component: " << vals(i,j,k)
00511                << " but reference D2 component: " << basisD2[l] << "\n";
00512             }
00513          }
00514       }
00515     }
00516 
00517     
00518     // Check D3 of basis function
00519     int D3cardinality = Intrepid::getDkCardinality(OPERATOR_D3, spaceDim);
00520     vals.resize(numFields, numPoints, D3cardinality);    
00521     hexBasis.getValues(vals, hexNodes, OPERATOR_D3);
00522      
00523     for (int i = 0; i < numFields; i++) {
00524       for (int j = 0; j < numPoints; j++) {
00525         for (int k = 0; k < D3cardinality; k++) {
00526           
00527           // basisD3 is (F,P,Dk), compute offset:
00528           int l = k + j * D3cardinality + i * D3cardinality * numPoints;
00529           if (std::abs(vals(i,j,k) - basisD3[l]) > INTREPID_TOL) {
00530             errorFlag++;
00531             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00532             
00533             // Output the multi-index of the value where the error is:
00534             *outStream << " At multi-index { ";
00535             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00536             *outStream << "}  computed D3 component: " << vals(i,j,k)
00537               << " but reference D3 component: " << basisD3[l] << "\n";
00538           }
00539         }
00540       }
00541     }
00542 
00543     
00544     // Check D4 of basis function
00545     int D4cardinality = Intrepid::getDkCardinality(OPERATOR_D4, spaceDim);
00546     vals.resize(numFields, numPoints, D4cardinality);    
00547     hexBasis.getValues(vals, hexNodes, OPERATOR_D4);
00548     for (int i = 0; i < numFields; i++) {
00549       for (int j = 0; j < numPoints; j++) {
00550         for (int k = 0; k < D4cardinality; k++) {
00551           
00552           // basisD4 is (F,P,Dk), compute offset:
00553           int l = k + j * D4cardinality + i * D4cardinality * numPoints;
00554           if (std::abs(vals(i,j,k) - basisD4[l]) > INTREPID_TOL) {
00555             errorFlag++;
00556             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00557             
00558             // Output the multi-index of the value where the error is:
00559             *outStream << " At multi-index { ";
00560             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00561             *outStream << "}  computed D4 component: " << vals(i,j,k)
00562               << " but reference D4 component: " << basisD2[l] << "\n";
00563           }
00564         }
00565       }
00566     }
00567     
00568 
00569 
00570   // Check D7 to D10 - must be zero. This basis does not cover D5 and D6
00571     for(EOperator op = OPERATOR_D7; op < OPERATOR_MAX; op++) {
00572       
00573       // The last dimension is the number of kth derivatives and needs to be resized for every Dk
00574       int DkCardin  = Intrepid::getDkCardinality(op, spaceDim);
00575       vals.resize(numFields, numPoints, DkCardin);    
00576 
00577       hexBasis.getValues(vals, hexNodes, op);
00578       for (int i = 0; i < vals.size(); i++) {
00579         if (std::abs(vals[i]) > INTREPID_TOL) {
00580           errorFlag++;
00581           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00582           
00583           // Get the multi-index of the value where the error is and the operator order
00584           std::vector<int> myIndex;
00585           vals.getMultiIndex(myIndex,i);
00586           int ord = Intrepid::getOperatorOrder(op);
00587           *outStream << " At multi-index { ";
00588           for(int j = 0; j < vals.rank(); j++) {
00589             *outStream << myIndex[j] << " ";
00590           }
00591           *outStream << "}  computed D"<< ord <<" component: " << vals[i] 
00592             << " but reference D" << ord << " component:  0 \n";
00593         }
00594       }
00595     }    
00596   }
00597   // Catch unexpected errors
00598   catch (std::logic_error err) {
00599     *outStream << err.what() << "\n\n";
00600     errorFlag = -1000;
00601   };
00602   
00603   if (errorFlag != 0)
00604     std::cout << "End Result: TEST FAILED\n";
00605   else
00606     std::cout << "End Result: TEST PASSED\n";
00607   
00608   // reset format state of std::cout
00609   std::cout.copyfmt(oldFormatState);
00610   
00611   return errorFlag;
00612 }