Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Basis/HGRAD_LINE_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 // 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 
00036 #include "Intrepid_FieldContainer.hpp"
00037 #include "Intrepid_HGRAD_LINE_C1_FEM.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_LINE_C1_FEM)                         |\n" \
00081     << "|                                                                             |\n" \
00082     << "|     1) Conversion of Dof tags into Dof ordinals and back                    |\n" \
00083     << "|     2) Basis values for VALUE, GRAD, CURL, DIV, and Dk operators            |\n" \
00084     << "|                                                                             |\n" \
00085     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00086     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00087     << "|                      Kara Peterson (dridzal@sandia.gov).                    |\n" \
00088     << "|                                                                             |\n" \
00089     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00090     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00091     << "|                                                                             |\n" \
00092     << "===============================================================================\n"\
00093     << "| TEST 1: Basis creation, exception testing                                   |\n"\
00094     << "===============================================================================\n";
00095 
00096   
00097   // Define basis and error flag
00098   Basis_HGRAD_LINE_C1_FEM<double, FieldContainer<double> > lineBasis;
00099   int errorFlag = 0;
00100 
00101   // Initialize throw counter for exception testing
00102   int nException     = 0;
00103   int throwCounter   = 0;
00104   
00105   // Define array containing the 2 vertices of the reference Line, its center and another point  
00106   FieldContainer<double> lineNodes(4, 1);
00107   lineNodes(0,0) =  -1.0;
00108   lineNodes(1,0) =   1.0;
00109   lineNodes(2,0) =   0.0;
00110   lineNodes(3,0) =   0.5;
00111 
00112   // Generic array for the output values; needs to be properly resized depending on the operator type
00113   FieldContainer<double> vals;
00114 
00115   try{
00116     // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
00117     vals.resize(lineBasis.getCardinality(), lineNodes.dimension(0) );
00118     
00119     // Exceptions 1-5: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 
00120     // getDofTag() to access invalid array elements thereby causing bounds check exception
00121     // exception #1
00122     INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(2,0,0), throwCounter, nException );
00123     // exception #2
00124     INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(1,1,1), throwCounter, nException );
00125     // exception #3
00126     INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(0,4,0), throwCounter, nException );
00127     // exception #4
00128     INTREPID_TEST_COMMAND( lineBasis.getDofTag(5), throwCounter, nException );
00129     // exception #5
00130     INTREPID_TEST_COMMAND( lineBasis.getDofTag(-1), throwCounter, nException );
00131 
00132 #ifdef HAVE_INTREPID_DEBUG
00133     // Exceptions 6-17 test exception handling with incorrectly dimensioned input/output arrays
00134     // exception #6: input points array must be of rank-2
00135     FieldContainer<double> badPoints1(4, 5, 3);
00136     INTREPID_TEST_COMMAND( lineBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00137     
00138     // exception #7 dimension 1 in the input point array must equal space dimension of the cell
00139     FieldContainer<double> badPoints2(4, 2);
00140     INTREPID_TEST_COMMAND( lineBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00141     
00142     // exception #8 output values must be of rank-2 for OPERATOR_VALUE
00143     FieldContainer<double> badVals1(4, 3, 1);
00144     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals1, lineNodes, OPERATOR_VALUE), throwCounter, nException );
00145     
00146     // exception #9 output values must be of rank-3 for OPERATOR_GRAD
00147     FieldContainer<double> badVals2(4, 3);
00148     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_GRAD), throwCounter, nException );
00149 
00150     // exception #10 output values must be of rank-3 for OPERATOR_DIV
00151     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_DIV), throwCounter, nException );
00152     
00153     // exception #11 output values must be of rank-3 for OPERATOR_CURL
00154     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_CURL), throwCounter, nException );
00155     
00156     // exception #12 output values must be of rank-3 for OPERATOR_D2
00157     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_D2), throwCounter, nException );
00158     
00159     // exception #13 incorrect 1st dimension of output array (must equal number of basis functions)
00160     FieldContainer<double> badVals3(lineBasis.getCardinality() + 1, lineNodes.dimension(0));
00161     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals3, lineNodes, OPERATOR_VALUE), throwCounter, nException );
00162     
00163     // exception #14 incorrect 0th dimension of output array (must equal number of points)
00164     FieldContainer<double> badVals4(lineBasis.getCardinality(), lineNodes.dimension(0) + 1);
00165     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals4, lineNodes, OPERATOR_VALUE), throwCounter, nException );
00166     
00167     // exception #15: incorrect 2nd dimension of output array (must equal the space dimension)
00168     FieldContainer<double> badVals5(lineBasis.getCardinality(), lineNodes.dimension(0), 2);
00169     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals5, lineNodes, OPERATOR_GRAD), throwCounter, nException );
00170     
00171     // exception #16: incorrect 2nd dimension of output array (must equal D2 cardinality in 1D)
00172     FieldContainer<double> badVals6(lineBasis.getCardinality(), lineNodes.dimension(0), 40);
00173     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals6, lineNodes, OPERATOR_D2), throwCounter, nException );
00174     
00175     // exception #17: incorrect 2nd dimension of output array (must equal D3 cardinality in 1D)
00176     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals6, lineNodes, OPERATOR_D3), throwCounter, nException );
00177 #endif
00178     
00179   }
00180   catch (std::logic_error err) {
00181     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00182     *outStream << err.what() << '\n';
00183     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00184     errorFlag = -1000;
00185   };
00186   
00187   // Check if number of thrown exceptions matches the one we expect
00188   if (throwCounter != nException) {
00189     errorFlag++;
00190     *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00191   }
00192   
00193   *outStream \
00194     << "\n"
00195     << "===============================================================================\n"\
00196     << "| TEST 2: correctness of tag to enum and enum to tag lookups                  |\n"\
00197     << "===============================================================================\n";
00198 
00199   try{
00200     std::vector<std::vector<int> > allTags = lineBasis.getAllDofTags();
00201     
00202     // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
00203     for (unsigned i = 0; i < allTags.size(); i++) {
00204       int bfOrd  = lineBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
00205       
00206       std::vector<int> myTag = lineBasis.getDofTag(bfOrd);
00207       if( !( (myTag[0] == allTags[i][0]) &&
00208              (myTag[1] == allTags[i][1]) &&
00209              (myTag[2] == allTags[i][2]) &&
00210              (myTag[3] == allTags[i][3]) ) ) {
00211         errorFlag++;
00212         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00213         *outStream << " getDofOrdinal( {" 
00214           << allTags[i][0] << ", " 
00215           << allTags[i][1] << ", " 
00216           << allTags[i][2] << ", " 
00217           << allTags[i][3] << "}) = " << bfOrd <<" but \n";   
00218         *outStream << " getDofTag(" << bfOrd << ") = { "
00219           << myTag[0] << ", " 
00220           << myTag[1] << ", " 
00221           << myTag[2] << ", " 
00222           << myTag[3] << "}\n";        
00223       }
00224     }
00225     
00226     // Now do the same but loop over basis functions
00227     for( int bfOrd = 0; bfOrd < lineBasis.getCardinality(); bfOrd++) {
00228       std::vector<int> myTag  = lineBasis.getDofTag(bfOrd);
00229       int myBfOrd = lineBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
00230       if( bfOrd != myBfOrd) {
00231         errorFlag++;
00232         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00233         *outStream << " getDofTag(" << bfOrd << ") = { "
00234           << myTag[0] << ", " 
00235           << myTag[1] << ", " 
00236           << myTag[2] << ", " 
00237           << myTag[3] << "} but getDofOrdinal({" 
00238           << myTag[0] << ", " 
00239           << myTag[1] << ", " 
00240           << myTag[2] << ", " 
00241           << myTag[3] << "} ) = " << myBfOrd << "\n";
00242       }
00243     }
00244   }
00245   catch (std::logic_error err){
00246     *outStream << err.what() << "\n\n";
00247     errorFlag = -1000;
00248   };
00249   
00250   *outStream \
00251     << "\n"
00252     << "===============================================================================\n"\
00253     << "| TEST 3: correctness of basis function values                                |\n"\
00254     << "===============================================================================\n";
00255 
00256   outStream -> precision(20);
00257 
00258   // VALUE: Each row gives the 2 correct basis set values at an evaluation point
00259   double basisValues[] = {
00260     1.0, 0.0,
00261     0.0, 1.0,
00262     0.5, 0.5,
00263     0.25, 0.75
00264   };
00265 
00266   // GRAD, DIV, CURL and D1: each row gives the 2 correct values of the gradients of the 2 basis functions
00267   double basisDerivs[] = {
00268     -0.5, 0.5,
00269     -0.5, 0.5,
00270     -0.5, 0.5,
00271     -0.5, 0.5
00272   };
00273 
00274   try{
00275     
00276     // Dimensions for the output arrays:
00277     int numFields = lineBasis.getCardinality();
00278     int numPoints = lineNodes.dimension(0);
00279     int spaceDim  = lineBasis.getBaseCellTopology().getDimension();
00280     
00281     // Generic array for values, grads, curls, etc. that will be properly sized before each call
00282     FieldContainer<double> vals;
00283     
00284     // Check VALUE of basis functions: resize vals to rank-2 container:
00285     vals.resize(numFields, numPoints);
00286     lineBasis.getValues(vals, lineNodes, OPERATOR_VALUE);
00287     for (int i = 0; i < numFields; i++) {
00288       for (int j = 0; j < numPoints; j++) {
00289           int l =  i + j * numFields;
00290            if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
00291              errorFlag++;
00292              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00293 
00294              // Output the multi-index of the value where the error is:
00295              *outStream << " At multi-index { ";
00296              *outStream << i << " ";*outStream << j << " ";
00297              *outStream << "}  computed value: " << vals(i,j)
00298                << " but reference value: " << basisValues[l] << "\n";
00299          }
00300       }
00301     }
00302     
00303     // Check derivatives of basis function: resize vals to rank-3 container
00304     vals.resize(numFields, numPoints, spaceDim);
00305     lineBasis.getValues(vals, lineNodes, OPERATOR_GRAD);
00306     for (int i = 0; i < numFields; i++) {
00307       for (int j = 0; j < numPoints; j++) {
00308         for (int k = 0; k < spaceDim; k++) {
00309            int l = k + i * spaceDim + j * spaceDim * numFields;
00310            if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
00311              errorFlag++;
00312              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00313 
00314              // Output the multi-index of the value where the error is:
00315              *outStream << " At multi-index { ";
00316              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00317              *outStream << "}  computed grad component: " << vals(i,j,k)
00318                << " but reference grad component: " << basisDerivs[l] << "\n";
00319             }
00320          }
00321       }
00322     }
00323 
00324     // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
00325     lineBasis.getValues(vals, lineNodes, OPERATOR_D1);
00326     for (int i = 0; i < numFields; i++) {
00327       for (int j = 0; j < numPoints; j++) {
00328         for (int k = 0; k < spaceDim; k++) {
00329            int l = k + i * spaceDim + j * spaceDim * numFields;
00330            if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
00331              errorFlag++;
00332              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00333 
00334              // Output the multi-index of the value where the error is:
00335              *outStream << " At multi-index { ";
00336              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00337              *outStream << "}  computed D1 component: " << vals(i,j,k)
00338                << " but reference D1 component: " << basisDerivs[l] << "\n";
00339             }
00340          }
00341       }
00342     }
00343 
00344     
00345     // Check CURL of basis function: resize vals just for illustration! 
00346     vals.resize(numFields, numPoints, spaceDim);
00347     lineBasis.getValues(vals, lineNodes, OPERATOR_CURL);
00348     for (int i = 0; i < numFields; i++) {
00349       for (int j = 0; j < numPoints; j++) {
00350         for (int k = 0; k < spaceDim; k++) {
00351            int l = k + i * spaceDim + j * spaceDim * numFields;
00352            if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
00353              errorFlag++;
00354              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00355 
00356              // Output the multi-index of the value where the error is:
00357              *outStream << " At multi-index { ";
00358              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00359              *outStream << "}  computed curl component: " << vals(i,j,k)
00360                << " but reference curl component: " << basisDerivs[l] << "\n";
00361             }
00362          }
00363       }
00364     }
00365 
00366 
00367     // Check DIV of basis function (do not resize vals because it has the correct size: DIV = CURL)
00368     lineBasis.getValues(vals, lineNodes, OPERATOR_DIV);
00369     for (int i = 0; i < numFields; i++) {
00370       for (int j = 0; j < numPoints; j++) {
00371         for (int k = 0; k < spaceDim; k++) {
00372            int l = k + i * spaceDim + j * spaceDim * numFields;
00373            if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
00374              errorFlag++;
00375              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00376 
00377              // Output the multi-index of the value where the error is:
00378              *outStream << " At multi-index { ";
00379              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00380              *outStream << "}  computed D1 component: " << vals(i,j,k)
00381                << " but reference D1 component: " << basisDerivs[l] << "\n";
00382             }
00383          }
00384       }
00385     }
00386 
00387 
00388     // Check all higher derivatives - must be zero. 
00389     for(EOperator op = OPERATOR_D2; op < OPERATOR_MAX; op++) {
00390       
00391       // The last dimension is the number of kth derivatives and needs to be resized for every Dk
00392       int DkCardin  = Intrepid::getDkCardinality(op, spaceDim);
00393       vals.resize(numFields, numPoints, DkCardin);    
00394       
00395       lineBasis.getValues(vals, lineNodes, op);
00396       for (int i = 0; i < vals.size(); i++) {
00397         if (std::abs(vals[i]) > INTREPID_TOL) {
00398           errorFlag++;
00399           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00400           
00401           // Get the multi-index of the value where the error is and the operator order
00402           std::vector<int> myIndex;
00403           vals.getMultiIndex(myIndex,i);
00404           int ord = Intrepid::getOperatorOrder(op);
00405           *outStream << " At multi-index { ";
00406           for(int j = 0; j < vals.rank(); j++) {
00407             *outStream << myIndex[j] << " ";
00408           }
00409           *outStream << "}  computed D"<< ord <<" component: " << vals[i] 
00410             << " but reference D" << ord << " component:  0 \n";
00411         }
00412       }
00413     }    
00414   }
00415 
00416   // Catch unexpected errors
00417   catch (std::logic_error err) {
00418     *outStream << err.what() << "\n\n";
00419     errorFlag = -1000;
00420   };
00421 
00422   if (errorFlag != 0)
00423     std::cout << "End Result: TEST FAILED\n";
00424   else
00425     std::cout << "End Result: TEST PASSED\n";
00426 
00427   // reset format state of std::cout
00428   std::cout.copyfmt(oldFormatState);
00429 
00430   return errorFlag;
00431 }