Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Basis/HGRAD_HEX_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_HEX_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_HEX_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_HEX_C2_FEM<double, FieldContainer<double> > hexBasis;
00097   int errorFlag = 0;
00098 
00099   // Initialize throw counter for exception testing
00100   int nException     = 0;
00101   int throwCounter   = 0;
00102 
00103   // Define arrayS containing the 27 nodes of hexahedron<27> topology
00104   FieldContainer<double> hexNodes(27, 3);
00105   // vertices
00106   hexNodes(0, 0) = -1.0;  hexNodes(0, 1) = -1.0;  hexNodes(0, 2) = -1.0;
00107   hexNodes(1, 0) =  1.0;  hexNodes(1, 1) = -1.0;  hexNodes(1, 2) = -1.0;
00108   hexNodes(2, 0) =  1.0;  hexNodes(2, 1) =  1.0;  hexNodes(2, 2) = -1.0;
00109   hexNodes(3, 0) = -1.0;  hexNodes(3, 1) =  1.0;  hexNodes(3, 2) = -1.0;
00110   
00111   hexNodes(4, 0) = -1.0;  hexNodes(4, 1) = -1.0;  hexNodes(4, 2) =  1.0;
00112   hexNodes(5, 0) =  1.0;  hexNodes(5, 1) = -1.0;  hexNodes(5, 2) =  1.0;
00113   hexNodes(6, 0) =  1.0;  hexNodes(6, 1) =  1.0;  hexNodes(6, 2) =  1.0;
00114   hexNodes(7, 0) = -1.0;  hexNodes(7, 1) =  1.0;  hexNodes(7, 2) =  1.0;  
00115   
00116   // nodes on edges
00117   hexNodes(8, 0) =  0.0;   hexNodes(8, 1) = -1.0;  hexNodes(8, 2) = -1.0;
00118   hexNodes(9, 0) =  1.0;   hexNodes(9, 1) =  0.0;  hexNodes(9, 2) = -1.0;
00119   hexNodes(10,0) =  0.0;   hexNodes(10,1) =  1.0;  hexNodes(10,2) = -1.0;
00120   hexNodes(11,0) = -1.0;   hexNodes(11,1) =  0.0;  hexNodes(11,2) = -1.0;
00121   hexNodes(12,0) = -1.0;   hexNodes(12,1) = -1.0;  hexNodes(12,2) =  0.0;
00122   hexNodes(13,0) =  1.0;   hexNodes(13,1) = -1.0;  hexNodes(13,2) =  0.0;
00123   hexNodes(14,0) =  1.0;   hexNodes(14,1) =  1.0;  hexNodes(14,2) =  0.0;
00124   hexNodes(15,0) = -1.0;   hexNodes(15,1) =  1.0;  hexNodes(15,2) =  0.0;
00125   hexNodes(16,0) =  0.0;   hexNodes(16,1) = -1.0;  hexNodes(16,2) =  1.0;
00126   hexNodes(17,0) =  1.0;   hexNodes(17,1) =  0.0;  hexNodes(17,2) =  1.0;
00127   hexNodes(18,0) =  0.0;   hexNodes(18,1) =  1.0;  hexNodes(18,2) =  1.0;
00128   hexNodes(19,0) = -1.0;   hexNodes(19,1) =  0.0;  hexNodes(19,2) =  1.0;
00129   
00130   // center
00131   hexNodes(20,0) =  0.0;  hexNodes(20,1) =  0.0;   hexNodes(20,2) =  0.0;
00132 
00133   // Face nodes
00134   hexNodes(21,0) =  0.0;   hexNodes(21,1) =  0.0;  hexNodes(21,2) = -1.0;
00135   hexNodes(22,0) =  0.0;   hexNodes(22,1) =  0.0;  hexNodes(22,2) =  1.0;
00136   hexNodes(23,0) = -1.0;   hexNodes(23,1) =  0.0;  hexNodes(23,2) =  0.0;
00137   hexNodes(24,0) =  1.0;   hexNodes(24,1) =  0.0;  hexNodes(24,2) =  0.0;
00138   hexNodes(25,0) =  0.0;   hexNodes(25,1) = -1.0;  hexNodes(25,2) =  0.0;
00139   hexNodes(26,0) =  0.0;   hexNodes(26,1) =  1.0;  hexNodes(26,2) =  0.0;
00140  
00141   // Generic array for the output values; needs to be properly resized depending on the operator type
00142   FieldContainer<double> vals;
00143 
00144   try{
00145     // exception #1: CURL cannot be applied to scalar functions in 3D
00146     // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
00147     vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0), 4 );
00148     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_CURL), throwCounter, nException );
00149 
00150     // exception #2: DIV cannot be applied to scalar functions in 3D
00151     // resize vals to rank-2 container with dimensions (num. basis functions, num. points)
00152     vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0) );
00153     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_DIV), throwCounter, nException );
00154         
00155     // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 
00156     // getDofTag() to access invalid array elements thereby causing bounds check exception
00157     // exception #3
00158     INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(3,10,0), throwCounter, nException );    
00159     // exception #4
00160     INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(1,2,1), throwCounter, nException );
00161     // exception #5
00162     INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(0,4,1), throwCounter, nException );
00163     // exception #6
00164     INTREPID_TEST_COMMAND( hexBasis.getDofTag(28), throwCounter, nException );
00165     // exception #7
00166     INTREPID_TEST_COMMAND( hexBasis.getDofTag(-1), throwCounter, nException );
00167 
00168 #ifdef HAVE_INTREPID_DEBUG 
00169     // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays
00170     // exception #8: input points array must be of rank-2
00171     FieldContainer<double> badPoints1(4, 5, 3);
00172     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00173 
00174     // exception #9 dimension 1 in the input point array must equal space dimension of the cell
00175     FieldContainer<double> badPoints2(4, hexBasis.getBaseCellTopology().getDimension() - 1);
00176     INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00177 
00178     // exception #10 output values must be of rank-2 for OPERATOR_VALUE
00179     FieldContainer<double> badVals1(4, 3, 1);
00180     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_VALUE), throwCounter, nException );
00181  
00182     // exception #11 output values must be of rank-3 for OPERATOR_GRAD
00183     FieldContainer<double> badVals2(4, 3);
00184     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_GRAD), throwCounter, nException );
00185     
00186     // exception #12 output values must be of rank-3 for OPERATOR_D1
00187     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D1), throwCounter, nException );
00188     
00189     // exception #13 output values must be of rank-3 for OPERATOR_D2
00190     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D2), throwCounter, nException );
00191     
00192     // exception #14 incorrect 0th dimension of output array (must equal number of basis functions)
00193     FieldContainer<double> badVals3(hexBasis.getCardinality() + 1, hexNodes.dimension(0));
00194     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals3, hexNodes, OPERATOR_VALUE), throwCounter, nException );
00195     
00196     // exception #15 incorrect 1st dimension of output array (must equal number of points)
00197     FieldContainer<double> badVals4(hexBasis.getCardinality(), hexNodes.dimension(0) + 1);
00198     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_VALUE), throwCounter, nException );
00199     
00200     // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
00201     FieldContainer<double> badVals5(hexBasis.getCardinality(), hexNodes.dimension(0), hexBasis.getBaseCellTopology().getDimension() - 1);
00202     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals5, hexNodes, OPERATOR_GRAD), throwCounter, nException );
00203     
00204     // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 3D)
00205     FieldContainer<double> badVals6(hexBasis.getCardinality(), hexNodes.dimension(0), 40);
00206     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals6, hexNodes, OPERATOR_D2), throwCounter, nException );
00207     
00208     // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 3D)
00209     FieldContainer<double> badVals7(hexBasis.getCardinality(), hexNodes.dimension(0), 50);
00210     INTREPID_TEST_COMMAND( hexBasis.getValues(badVals7, hexNodes, OPERATOR_D3), throwCounter, nException );
00211 #endif
00212 
00213   }
00214   catch (std::logic_error err) {
00215     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00216     *outStream << err.what() << '\n';
00217     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00218     errorFlag = -1000;
00219   };
00220   
00221   // Check if number of thrown exceptions matches the one we expect 
00222   // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match.
00223   if (throwCounter != nException) {
00224     errorFlag++;
00225     *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00226   }
00227   
00228   *outStream \
00229     << "\n"
00230     << "===============================================================================\n"\
00231     << "| TEST 2: correctness of tag to enum and enum to tag lookups                  |\n"\
00232     << "===============================================================================\n";
00233   
00234   try{
00235     std::vector<std::vector<int> > allTags = hexBasis.getAllDofTags();
00236     
00237     // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
00238     for (unsigned i = 0; i < allTags.size(); i++) {
00239       int bfOrd  = hexBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
00240       
00241       std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
00242        if( !( (myTag[0] == allTags[i][0]) &&
00243               (myTag[1] == allTags[i][1]) &&
00244               (myTag[2] == allTags[i][2]) &&
00245               (myTag[3] == allTags[i][3]) ) ) {
00246         errorFlag++;
00247         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00248         *outStream << " getDofOrdinal( {" 
00249           << allTags[i][0] << ", " 
00250           << allTags[i][1] << ", " 
00251           << allTags[i][2] << ", " 
00252           << allTags[i][3] << "}) = " << bfOrd <<" but \n";   
00253         *outStream << " getDofTag(" << bfOrd << ") = { "
00254           << myTag[0] << ", " 
00255           << myTag[1] << ", " 
00256           << myTag[2] << ", " 
00257           << myTag[3] << "}\n";        
00258       }
00259     }
00260     
00261     // Now do the same but loop over basis functions
00262     for( int bfOrd = 0; bfOrd < hexBasis.getCardinality(); bfOrd++) {
00263       std::vector<int> myTag  = hexBasis.getDofTag(bfOrd);
00264       int myBfOrd = hexBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
00265       if( bfOrd != myBfOrd) {
00266         errorFlag++;
00267         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00268         *outStream << " getDofTag(" << bfOrd << ") = { "
00269           << myTag[0] << ", " 
00270           << myTag[1] << ", " 
00271           << myTag[2] << ", " 
00272           << myTag[3] << "} but getDofOrdinal({" 
00273           << myTag[0] << ", " 
00274           << myTag[1] << ", " 
00275           << myTag[2] << ", " 
00276           << myTag[3] << "} ) = " << myBfOrd << "\n";
00277       }
00278     }
00279   }
00280   catch (std::logic_error err){
00281     *outStream << err.what() << "\n\n";
00282     errorFlag = -1000;
00283   };
00284   
00285   
00286   *outStream \
00287     << "\n"
00288     << "===============================================================================\n"\
00289     << "| TEST 3: correctness of basis function values                                |\n"\
00290     << "===============================================================================\n";
00291   
00292   outStream -> precision(20);
00293   
00294   // VALUE: Each row gives the 8 correct basis set values at an evaluation point
00295   double basisValues[] = {
00296     1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00297     0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00298     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00299     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, \
00300     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00301     0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00302     0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00303     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, \
00304     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, \
00305     0, 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, 1.000, 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, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00308     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, \
00309     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, \
00310     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00311     0, 0, 0, 1.000, 0, 0, 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, 1.000, 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, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, \
00314     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, \
00315     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00316     0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00317     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 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, 1.000, 0, 0, 0, 0, 0, \
00319     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00320     1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00321     0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00322     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00323     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 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, \
00325     0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00326     0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
00327     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000, 0, 0, 0, 0, 0, 0, 0, 0, \
00328     0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.000  };
00329   
00330   
00331   // GRAD, D1, D2, D3 and D4 test values are stored in files due to their large size
00332   std::string     fileName;
00333   std::ifstream   dataFile;
00334 
00335   // GRAD and D1 values are stored in (F,P,D) format in a data file. Read file and do the test
00336   std::vector<double> basisGrads;           // Flat array for the gradient values.
00337   
00338   fileName = "./testdata/HEX_C2_GradVals.dat";
00339   dataFile.open(fileName.c_str());
00340   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00341                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open GRAD values data file, test aborted.");
00342   while (!dataFile.eof() ){
00343     double temp;
00344     string line;                            // string for one line of input file
00345     std::getline(dataFile, line);           // get next line from file
00346     stringstream data_line(line);           // convert to stringstream
00347     while(data_line >> temp){               // extract value from line
00348       basisGrads.push_back(temp);           // push into vector
00349     }
00350   }
00351   // It turns out that just closing and then opening the ifstream variable does not reset it
00352   // and subsequent open() command fails. One fix is to explicitely clear the ifstream, or
00353   // scope the variables.
00354   dataFile.close();
00355   dataFile.clear();
00356     
00357    
00358   //D2: flat array with the values of D2 applied to basis functions. Multi-index is (F,P,D2cardinality)
00359   std::vector<double> basisD2; 
00360   fileName = "./testdata/HEX_C2_D2Vals.dat";  
00361   dataFile.open(fileName.c_str());
00362   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00363                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open D2 values data file, test aborted.");
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       basisD2.push_back(temp);              // push into vector
00371     }
00372   }
00373   dataFile.close();
00374   dataFile.clear();
00375 
00376   
00377   //D3: flat array with the values of D3 applied to basis functions. Multi-index is (F,P,D3cardinality)
00378   std::vector<double> basisD3;
00379   
00380   fileName = "./testdata/HEX_C2_D3Vals.dat";  
00381   dataFile.open(fileName.c_str());
00382   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00383                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open D3 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       basisD3.push_back(temp);              // push into vector
00392     }
00393   }
00394   dataFile.close();
00395   dataFile.clear();
00396  
00397   
00398   //D4: flat array with the values of D3 applied to basis functions. Multi-index is (F,P,D4cardinality)
00399   std::vector<double> basisD4;
00400   
00401   fileName = "./testdata/HEX_C2_D4Vals.dat";
00402   dataFile.open(fileName.c_str());
00403   TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
00404                       ">>> ERROR (HGRAD_HEX_C2/test01): could not open D4 values data file, test aborted.");
00405   
00406   while (!dataFile.eof() ){
00407     double temp;
00408     string line;                            // string for one line of input file
00409     std::getline(dataFile, line);           // get next line from file
00410     stringstream data_line(line);           // convert to stringstream
00411     while(data_line >> temp){               // extract value from line
00412       basisD4.push_back(temp);              // push into vector
00413     }
00414   }
00415   dataFile.close();
00416   dataFile.clear();
00417   
00418 
00419   try{
00420         
00421     // Dimensions for the output arrays:
00422     int numFields = hexBasis.getCardinality();
00423     int numPoints = hexNodes.dimension(0);
00424     int spaceDim  = hexBasis.getBaseCellTopology().getDimension();
00425     
00426     // Generic array for values, grads, curls, etc. that will be properly sized before each call
00427     FieldContainer<double> vals;
00428     
00429     // Check VALUE of basis functions: resize vals to rank-2 container:
00430     vals.resize(numFields, numPoints);
00431     hexBasis.getValues(vals, hexNodes, OPERATOR_VALUE);
00432     for (int i = 0; i < numFields; i++) {
00433       for (int j = 0; j < numPoints; j++) {
00434           int l =  i + j * numFields;
00435            if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
00436              errorFlag++;
00437              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00438 
00439              // Output the multi-index of the value where the error is:
00440              *outStream << " At multi-index { ";
00441              *outStream << i << " ";*outStream << j << " ";
00442              *outStream << "}  computed value: " << vals(i,j)
00443                << " but reference value: " << basisValues[l] << "\n";
00444          }
00445       }
00446     }
00447     
00448     
00449     // Check GRAD of basis function: resize vals to rank-3 container
00450     vals.resize(numFields, numPoints, spaceDim);
00451     hexBasis.getValues(vals, hexNodes, OPERATOR_GRAD);
00452     for (int i = 0; i < numFields; i++) {
00453       for (int j = 0; j < numPoints; j++) {
00454         for (int k = 0; k < spaceDim; k++) {
00455           
00456           // basisGrads is (F,P,D), compute offset:
00457           int l = k + j * spaceDim + i * spaceDim * numPoints;
00458            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00459              errorFlag++;
00460              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00461 
00462              // Output the multi-index of the value where the error is:
00463              *outStream << " At multi-index { ";
00464              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00465              *outStream << "}  computed grad component: " << vals(i,j,k)
00466                << " but reference grad component: " << basisGrads[l] << "\n";
00467             }
00468          }
00469       }
00470     }
00471     
00472     // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
00473     hexBasis.getValues(vals, hexNodes, OPERATOR_D1);
00474     for (int i = 0; i < numFields; i++) {
00475       for (int j = 0; j < numPoints; j++) {
00476         for (int k = 0; k < spaceDim; k++) {
00477 
00478           // basisGrads is (F,P,D), compute offset:
00479           int l = k + j * spaceDim + i * spaceDim * numPoints;
00480            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00481              errorFlag++;
00482              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00483 
00484              // Output the multi-index of the value where the error is:
00485              *outStream << " At multi-index { ";
00486              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00487              *outStream << "}  computed D1 component: " << vals(i,j,k)
00488                << " but reference D1 component: " << basisGrads[l] << "\n";
00489             }
00490          }
00491       }
00492     }
00493 
00494     
00495     // Check D2 of basis function
00496     int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim);
00497     vals.resize(numFields, numPoints, D2cardinality);    
00498     hexBasis.getValues(vals, hexNodes, OPERATOR_D2);
00499     for (int i = 0; i < numFields; i++) {
00500       for (int j = 0; j < numPoints; j++) {
00501         for (int k = 0; k < D2cardinality; k++) {
00502  
00503           // basisD2 is (F,P,Dk), compute offset:
00504           int l = k + j * D2cardinality + i * D2cardinality * numPoints;
00505            if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) {
00506              errorFlag++;
00507              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00508 
00509              // Output the multi-index of the value where the error is:
00510              *outStream << " At multi-index { ";
00511              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00512              *outStream << "}  computed D2 component: " << vals(i,j,k)
00513                << " but reference D2 component: " << basisD2[l] << "\n";
00514             }
00515          }
00516       }
00517     }
00518 
00519     
00520     // Check D3 of basis function
00521     int D3cardinality = Intrepid::getDkCardinality(OPERATOR_D3, spaceDim);
00522     vals.resize(numFields, numPoints, D3cardinality);    
00523     hexBasis.getValues(vals, hexNodes, OPERATOR_D3);
00524      
00525     for (int i = 0; i < numFields; i++) {
00526       for (int j = 0; j < numPoints; j++) {
00527         for (int k = 0; k < D3cardinality; k++) {
00528           
00529           // basisD3 is (F,P,Dk), compute offset:
00530           int l = k + j * D3cardinality + i * D3cardinality * numPoints;
00531           if (std::abs(vals(i,j,k) - basisD3[l]) > INTREPID_TOL) {
00532             errorFlag++;
00533             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00534             
00535             // Output the multi-index of the value where the error is:
00536             *outStream << " At multi-index { ";
00537             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00538             *outStream << "}  computed D3 component: " << vals(i,j,k)
00539               << " but reference D3 component: " << basisD3[l] << "\n";
00540           }
00541         }
00542       }
00543     }
00544 
00545     
00546     // Check D4 of basis function
00547     int D4cardinality = Intrepid::getDkCardinality(OPERATOR_D4, spaceDim);
00548     vals.resize(numFields, numPoints, D4cardinality);    
00549     hexBasis.getValues(vals, hexNodes, OPERATOR_D4);
00550     for (int i = 0; i < numFields; i++) {
00551       for (int j = 0; j < numPoints; j++) {
00552         for (int k = 0; k < D4cardinality; k++) {
00553           
00554           // basisD4 is (F,P,Dk), compute offset:
00555           int l = k + j * D4cardinality + i * D4cardinality * numPoints;
00556           if (std::abs(vals(i,j,k) - basisD4[l]) > INTREPID_TOL) {
00557             errorFlag++;
00558             *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00559             
00560             // Output the multi-index of the value where the error is:
00561             *outStream << " At multi-index { ";
00562             *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00563             *outStream << "}  computed D4 component: " << vals(i,j,k)
00564               << " but reference D4 component: " << basisD2[l] << "\n";
00565           }
00566         }
00567       }
00568     }
00569     
00570 
00571 
00572     // Check D7 to D10 - must be zero. This basis does not support D5 and D6
00573     for(EOperator op = OPERATOR_D7; op < OPERATOR_MAX; op++) {
00574       
00575       // The last dimension is the number of kth derivatives and needs to be resized for every Dk
00576       int DkCardin  = Intrepid::getDkCardinality(op, spaceDim);
00577       vals.resize(numFields, numPoints, DkCardin);    
00578 
00579       hexBasis.getValues(vals, hexNodes, op);
00580       for (int i = 0; i < vals.size(); i++) {
00581         if (std::abs(vals[i]) > INTREPID_TOL) {
00582           errorFlag++;
00583           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00584           
00585           // Get the multi-index of the value where the error is and the operator order
00586           std::vector<int> myIndex;
00587           vals.getMultiIndex(myIndex,i);
00588           int ord = Intrepid::getOperatorOrder(op);
00589           *outStream << " At multi-index { ";
00590           for(int j = 0; j < vals.rank(); j++) {
00591             *outStream << myIndex[j] << " ";
00592           }
00593           *outStream << "}  computed D"<< ord <<" component: " << vals[i] 
00594             << " but reference D" << ord << " component:  0 \n";
00595         }
00596       }
00597     }    
00598   }
00599   
00600   // Catch unexpected errors
00601   catch (std::logic_error err) {
00602     *outStream << err.what() << "\n\n";
00603     errorFlag = -1000;
00604   };
00605 
00606   *outStream \
00607     << "\n"
00608     << "===============================================================================\n"\
00609     << "| TEST 4: correctness of DoF locations                                        |\n"\
00610     << "===============================================================================\n";
00611 
00612   try{
00613     Teuchos::RCP<Basis<double, FieldContainer<double> > > basis =
00614       Teuchos::rcp(new Basis_HGRAD_HEX_C2_FEM<double, FieldContainer<double> >);
00615     Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface =
00616       Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis);
00617 
00618     FieldContainer<double> cvals;
00619     FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality());
00620 
00621     // Check exceptions.
00622 #ifdef HAVE_INTREPID_DEBUG
00623     cvals.resize(1,2,3);
00624     INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
00625     cvals.resize(3,2);
00626     INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
00627     cvals.resize(27,2);
00628     INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
00629 #endif
00630     cvals.resize(27,3);
00631     INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--;
00632     // Check if number of thrown exceptions matches the one we expect
00633     if (throwCounter != nException) {
00634       errorFlag++;
00635       *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00636     }
00637 
00638     // Check mathematical correctness.
00639     basis->getValues(bvals, cvals, OPERATOR_VALUE);
00640     char buffer[120];
00641     for (int i=0; i<bvals.dimension(0); i++) {
00642       for (int j=0; j<bvals.dimension(1); j++) {
00643         if ((i != j) && (std::abs(bvals(i,j) - 0.0) > INTREPID_TOL)) {
00644           errorFlag++;
00645           sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), cvals(i,2), bvals(i,j), 0.0);
00646           *outStream << buffer;
00647         }
00648         else if ((i == j) && (std::abs(bvals(i,j) - 1.0) > INTREPID_TOL)) {
00649           errorFlag++;
00650           sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), cvals(i,2), bvals(i,j), 1.0);
00651           *outStream << buffer;
00652         }
00653       }
00654     }
00655 
00656   }
00657   catch (std::logic_error err){
00658     *outStream << err.what() << "\n\n";
00659     errorFlag = -1000;
00660   };
00661   
00662   if (errorFlag != 0)
00663     std::cout << "End Result: TEST FAILED\n";
00664   else
00665     std::cout << "End Result: TEST PASSED\n";
00666   
00667   // reset format state of std::cout
00668   std::cout.copyfmt(oldFormatState);
00669   
00670   return errorFlag;
00671 }