Intrepid
http://trilinos.sandia.gov/packages/docs/r10.12/packages/intrepid/test/Discretization/Basis/HGRAD_TRI_C1_FEM/test_01.cpp
00001 // @HEADER
00002 // ************************************************************************
00003 //
00004 //                           Intrepid Package
00005 //                 Copyright (2007) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Pavel Bochev  (pbboche@sandia.gov)
00038 //                    Denis Ridzal  (dridzal@sandia.gov), or
00039 //                    Kara Peterson (kjpeter@sandia.gov)
00040 //
00041 // ************************************************************************
00042 // @HEADER
00043 
00049 #include "Intrepid_FieldContainer.hpp"
00050 #include "Intrepid_HGRAD_TRI_C1_FEM.hpp"
00051 #include "Teuchos_oblackholestream.hpp"
00052 #include "Teuchos_RCP.hpp"
00053 #include "Teuchos_GlobalMPISession.hpp"
00054 
00055 using namespace std;
00056 using namespace Intrepid;
00057 
00058 #define INTREPID_TEST_COMMAND( S , throwCounter, nException )                                                              \
00059 {                                                                                                                          \
00060   ++nException;                                                                                                            \
00061   try {                                                                                                                    \
00062     S ;                                                                                                                    \
00063   }                                                                                                                        \
00064   catch (std::logic_error err) {                                                                                           \
00065       ++throwCounter;                                                                                                      \
00066       *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
00067       *outStream << err.what() << '\n';                                                                                    \
00068       *outStream << "-------------------------------------------------------------------------------" << "\n\n";           \
00069   };                                                                                                                       \
00070 }
00071 int main(int argc, char *argv[]) {
00072 
00073   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00074 
00075   // This little trick lets us print to std::cout only if
00076   // a (dummy) command-line argument is provided.
00077   int iprint     = argc - 1;
00078   Teuchos::RCP<std::ostream> outStream;
00079   Teuchos::oblackholestream bhs; // outputs nothing
00080   if (iprint > 0)
00081     outStream = Teuchos::rcp(&std::cout, false);
00082   else
00083     outStream = Teuchos::rcp(&bhs, false);
00084 
00085   // Save the format state of the original std::cout.
00086   Teuchos::oblackholestream oldFormatState;
00087   oldFormatState.copyfmt(std::cout);
00088 
00089   *outStream \
00090     << "===============================================================================\n" \
00091     << "|                                                                             |\n" \
00092     << "|                 Unit Test (Basis_HGRAD_TRI_C1_FEM)                          |\n" \
00093     << "|                                                                             |\n" \
00094     << "|     1) Conversion of Dof tags into Dof ordinals and back                    |\n" \
00095     << "|     2) Basis values for VALUE, GRAD, CURL, and Dk operators                 |\n" \
00096     << "|                                                                             |\n" \
00097     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00098     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00099     << "|                      Kara Peterson (dridzal@sandia.gov).                    |\n" \
00100     << "|                                                                             |\n" \
00101     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00102     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00103     << "|                                                                             |\n" \
00104     << "===============================================================================\n"\
00105     << "| TEST 1: Basis creation, exception testing                                   |\n"\
00106     << "===============================================================================\n";
00107 
00108   
00109   // Define basis and error flag
00110   Basis_HGRAD_TRI_C1_FEM<double, FieldContainer<double> > triBasis;
00111   int errorFlag = 0;
00112 
00113   // Initialize throw counter for exception testing
00114   int nException     = 0;
00115   int throwCounter   = 0;
00116   
00117   // Define array containing the 3 vertices of the reference Triangle, its center and another point  
00118   FieldContainer<double> triNodes(5, 2);
00119   triNodes(0,0) =  0.0;  triNodes(0,1) =  0.0;
00120   triNodes(1,0) =  1.0;  triNodes(1,1) =  0.0;
00121   triNodes(2,0) =  0.0;  triNodes(2,1) =  1.0;
00122   triNodes(3,0) =  0.5;  triNodes(3,1) =  0.5;
00123   triNodes(4,0) =  0.0;  triNodes(4,1) =  0.75;
00124   
00125   // Generic array for the output values; needs to be properly resized depending on the operator type
00126   FieldContainer<double> vals;
00127 
00128   try{
00129     // exception #1: DIV cannot be applied to scalar functions
00130     // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
00131     vals.resize(triBasis.getCardinality(), triNodes.dimension(0) );
00132     INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_DIV), throwCounter, nException );
00133     
00134     // Exceptions 2-6: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 
00135     // getDofTag() to access invalid array elements thereby causing bounds check exception
00136     // exception #2
00137     INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(2,0,0), throwCounter, nException );
00138     // exception #3
00139     INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(1,1,1), throwCounter, nException );
00140     // exception #4
00141     INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(0,4,0), throwCounter, nException );
00142     // exception #5
00143     INTREPID_TEST_COMMAND( triBasis.getDofTag(5), throwCounter, nException );
00144     // exception #6
00145     INTREPID_TEST_COMMAND( triBasis.getDofTag(-1), throwCounter, nException );
00146 
00147 #ifdef HAVE_INTREPID_DEBUG
00148     // Exceptions 7-17 test exception handling with incorrectly dimensioned input/output arrays
00149     // exception #7: input points array must be of rank-2
00150     FieldContainer<double> badPoints1(4, 5, 3);
00151     INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00152     
00153     // exception #8 dimension 1 in the input point array must equal space dimension of the cell
00154     FieldContainer<double> badPoints2(4, 3);
00155     INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00156     
00157     // exception #9 output values must be of rank-2 for OPERATOR_VALUE
00158     FieldContainer<double> badVals1(4, 3, 1);
00159     INTREPID_TEST_COMMAND( triBasis.getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException );
00160     
00161     // exception #10 output values must be of rank-3 for OPERATOR_GRAD
00162     FieldContainer<double> badVals2(4, 3);
00163     INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_GRAD), throwCounter, nException );
00164     
00165     // exception #11 output values must be of rank-3 for OPERATOR_CURL
00166     INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_CURL), throwCounter, nException );
00167     
00168     // exception #12 output values must be of rank-3 for OPERATOR_D2
00169     INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_D2), throwCounter, nException );
00170     
00171     // exception #13 incorrect 1st dimension of output array (must equal number of basis functions)
00172     FieldContainer<double> badVals3(triBasis.getCardinality() + 1, triNodes.dimension(0));
00173     INTREPID_TEST_COMMAND( triBasis.getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException );
00174     
00175     // exception #14 incorrect 0th dimension of output array (must equal number of points)
00176     FieldContainer<double> badVals4(triBasis.getCardinality(), triNodes.dimension(0) + 1);
00177     INTREPID_TEST_COMMAND( triBasis.getValues(badVals4, triNodes, OPERATOR_VALUE), throwCounter, nException );
00178     
00179     // exception #15: incorrect 2nd dimension of output array (must equal the space dimension)
00180     FieldContainer<double> badVals5(triBasis.getCardinality(), triNodes.dimension(0), 4);
00181     INTREPID_TEST_COMMAND( triBasis.getValues(badVals5, triNodes, OPERATOR_GRAD), throwCounter, nException );
00182     
00183     // exception #16: incorrect 2nd dimension of output array (must equal D2 cardinality in 2D)
00184     FieldContainer<double> badVals6(triBasis.getCardinality(), triNodes.dimension(0), 40);
00185     INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_D2), throwCounter, nException );
00186     
00187     // exception #17: incorrect 2nd dimension of output array (must equal D3 cardinality in 2D)
00188     INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_D3), throwCounter, nException );
00189 #endif
00190     
00191   }
00192   catch (std::logic_error err) {
00193     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00194     *outStream << err.what() << '\n';
00195     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00196     errorFlag = -1000;
00197   };
00198   
00199   // Check if number of thrown exceptions matches the one we expect
00200   if (throwCounter != nException) {
00201     errorFlag++;
00202     *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00203   }
00204   
00205   *outStream \
00206     << "\n"
00207     << "===============================================================================\n"\
00208     << "| TEST 2: correctness of tag to enum and enum to tag lookups                  |\n"\
00209     << "===============================================================================\n";
00210 
00211   try{
00212     std::vector<std::vector<int> > allTags = triBasis.getAllDofTags();
00213     
00214     // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
00215     for (unsigned i = 0; i < allTags.size(); i++) {
00216       int bfOrd  = triBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
00217       
00218       std::vector<int> myTag = triBasis.getDofTag(bfOrd);
00219       if( !( (myTag[0] == allTags[i][0]) &&
00220              (myTag[1] == allTags[i][1]) &&
00221              (myTag[2] == allTags[i][2]) &&
00222              (myTag[3] == allTags[i][3]) ) ) {
00223         errorFlag++;
00224         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00225         *outStream << " getDofOrdinal( {" 
00226           << allTags[i][0] << ", " 
00227           << allTags[i][1] << ", " 
00228           << allTags[i][2] << ", " 
00229           << allTags[i][3] << "}) = " << bfOrd <<" but \n";   
00230         *outStream << " getDofTag(" << bfOrd << ") = { "
00231           << myTag[0] << ", " 
00232           << myTag[1] << ", " 
00233           << myTag[2] << ", " 
00234           << myTag[3] << "}\n";        
00235       }
00236     }
00237     
00238     // Now do the same but loop over basis functions
00239     for( int bfOrd = 0; bfOrd < triBasis.getCardinality(); bfOrd++) {
00240       std::vector<int> myTag  = triBasis.getDofTag(bfOrd);
00241       int myBfOrd = triBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
00242       if( bfOrd != myBfOrd) {
00243         errorFlag++;
00244         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00245         *outStream << " getDofTag(" << bfOrd << ") = { "
00246           << myTag[0] << ", " 
00247           << myTag[1] << ", " 
00248           << myTag[2] << ", " 
00249           << myTag[3] << "} but getDofOrdinal({" 
00250           << myTag[0] << ", " 
00251           << myTag[1] << ", " 
00252           << myTag[2] << ", " 
00253           << myTag[3] << "} ) = " << myBfOrd << "\n";
00254       }
00255     }
00256   }
00257   catch (std::logic_error err){
00258     *outStream << err.what() << "\n\n";
00259     errorFlag = -1000;
00260   };
00261   
00262   *outStream \
00263     << "\n"
00264     << "===============================================================================\n"\
00265     << "| TEST 3: correctness of basis function values                                |\n"\
00266     << "===============================================================================\n";
00267 
00268   outStream -> precision(20);
00269 
00270   // VALUE: Each row gives the 3 correct basis set values at an evaluation point
00271   double basisValues[] = {
00272     1.0, 0.0, 0.0,
00273     0.0, 1.0, 0.0,
00274     0.0, 0.0, 1.0,
00275     0.0, 0.5, 0.5,
00276     0.25,0.0, 0.75
00277   };
00278 
00279   // GRAD and D1: each row gives the 6 correct values of the gradients of the 3 basis functions
00280   double basisGrads[] = {
00281     -1.0, -1.0,    1.0,  0.0,    0.0,  1.0,
00282     -1.0, -1.0,    1.0,  0.0,    0.0,  1.0,
00283     -1.0, -1.0,    1.0,  0.0,    0.0,  1.0,
00284     -1.0, -1.0,    1.0,  0.0,    0.0,  1.0,
00285     -1.0, -1.0,    1.0,  0.0,    0.0,  1.0,
00286   };
00287 
00288   // CURL: each row gives the 6 correct values of the curls of the 3 basis functions
00289   double basisCurls[] = {
00290     -1.0,  1.0,    0.0, -1.0,    1.0,  0.0,
00291     -1.0,  1.0,    0.0, -1.0,    1.0,  0.0,
00292     -1.0,  1.0,    0.0, -1.0,    1.0,  0.0,
00293     -1.0,  1.0,    0.0, -1.0,    1.0,  0.0,
00294     -1.0,  1.0,    0.0, -1.0,    1.0,  0.0
00295   };
00296 
00297   try{
00298     
00299     // Dimensions for the output arrays:
00300     int numFields = triBasis.getCardinality();
00301     int numPoints = triNodes.dimension(0);
00302     int spaceDim  = triBasis.getBaseCellTopology().getDimension();
00303     
00304     // Generic array for values, grads, curls, etc. that will be properly sized before each call
00305     FieldContainer<double> vals;
00306     
00307     // Check VALUE of basis functions: resize vals to rank-2 container:
00308     vals.resize(numFields, numPoints);
00309     triBasis.getValues(vals, triNodes, OPERATOR_VALUE);
00310     for (int i = 0; i < numFields; i++) {
00311       for (int j = 0; j < numPoints; j++) {
00312           int l =  i + j * numFields;
00313            if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
00314              errorFlag++;
00315              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00316 
00317              // Output the multi-index of the value where the error is:
00318              *outStream << " At multi-index { ";
00319              *outStream << i << " ";*outStream << j << " ";
00320              *outStream << "}  computed value: " << vals(i,j)
00321                << " but reference value: " << basisValues[l] << "\n";
00322          }
00323       }
00324     }
00325     
00326     // Check GRAD of basis function: resize vals to rank-3 container
00327     vals.resize(numFields, numPoints, spaceDim);
00328     triBasis.getValues(vals, triNodes, OPERATOR_GRAD);
00329     for (int i = 0; i < numFields; i++) {
00330       for (int j = 0; j < numPoints; j++) {
00331         for (int k = 0; k < spaceDim; k++) {
00332            int l = k + i * spaceDim + j * spaceDim * numFields;
00333            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00334              errorFlag++;
00335              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00336 
00337              // Output the multi-index of the value where the error is:
00338              *outStream << " At multi-index { ";
00339              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00340              *outStream << "}  computed grad component: " << vals(i,j,k)
00341                << " but reference grad component: " << basisGrads[l] << "\n";
00342             }
00343          }
00344       }
00345     }
00346 
00347     // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
00348     triBasis.getValues(vals, triNodes, OPERATOR_D1);
00349     for (int i = 0; i < numFields; i++) {
00350       for (int j = 0; j < numPoints; j++) {
00351         for (int k = 0; k < spaceDim; k++) {
00352            int l = k + i * spaceDim + j * spaceDim * numFields;
00353            if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
00354              errorFlag++;
00355              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00356 
00357              // Output the multi-index of the value where the error is:
00358              *outStream << " At multi-index { ";
00359              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00360              *outStream << "}  computed D1 component: " << vals(i,j,k)
00361                << " but reference D1 component: " << basisGrads[l] << "\n";
00362             }
00363          }
00364       }
00365     }
00366 
00367     
00368     // Check CURL of basis function: resize vals just for illustration! 
00369     vals.resize(numFields, numPoints, spaceDim);
00370     triBasis.getValues(vals, triNodes, OPERATOR_CURL);
00371     for (int i = 0; i < numFields; i++) {
00372       for (int j = 0; j < numPoints; j++) {
00373         for (int k = 0; k < spaceDim; k++) {
00374            int l = k + i * spaceDim + j * spaceDim * numFields;
00375            if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) {
00376              errorFlag++;
00377              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00378 
00379              // Output the multi-index of the value where the error is:
00380              *outStream << " At multi-index { ";
00381              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00382              *outStream << "}  computed curl component: " << vals(i,j,k)
00383                << " but reference curl component: " << basisCurls[l] << "\n";
00384             }
00385          }
00386       }
00387     }
00388 
00389     // Check all higher derivatives - must be zero. 
00390     for(EOperator op = OPERATOR_D2; op < OPERATOR_MAX; op++) {
00391       
00392       // The last dimension is the number of kth derivatives and needs to be resized for every Dk
00393       int DkCardin  = Intrepid::getDkCardinality(op, spaceDim);
00394       vals.resize(numFields, numPoints, DkCardin);    
00395       
00396       triBasis.getValues(vals, triNodes, op);
00397       for (int i = 0; i < vals.size(); i++) {
00398         if (std::abs(vals[i]) > INTREPID_TOL) {
00399           errorFlag++;
00400           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00401           
00402           // Get the multi-index of the value where the error is and the operator order
00403           std::vector<int> myIndex;
00404           vals.getMultiIndex(myIndex,i);
00405           int ord = Intrepid::getOperatorOrder(op);
00406           *outStream << " At multi-index { ";
00407           for(int j = 0; j < vals.rank(); j++) {
00408             *outStream << myIndex[j] << " ";
00409           }
00410           *outStream << "}  computed D"<< ord <<" component: " << vals[i] 
00411             << " but reference D" << ord << " component:  0 \n";
00412         }
00413       }
00414     }    
00415   }
00416 
00417   // Catch unexpected errors
00418   catch (std::logic_error err) {
00419     *outStream << err.what() << "\n\n";
00420     errorFlag = -1000;
00421   };
00422 
00423   if (errorFlag != 0)
00424     std::cout << "End Result: TEST FAILED\n";
00425   else
00426     std::cout << "End Result: TEST PASSED\n";
00427 
00428   // reset format state of std::cout
00429   std::cout.copyfmt(oldFormatState);
00430 
00431   return errorFlag;
00432 }