Intrepid
http://trilinos.sandia.gov/packages/docs/r10.12/packages/intrepid/test/Discretization/Basis/HGRAD_LINE_Cn_FEM_JACOBI/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_LINE_Cn_FEM_JACOBI.hpp"
00051 #include "Intrepid_DefaultCubatureFactory.hpp"
00052 #include "Intrepid_ArrayTools.hpp"
00053 #include "Intrepid_FunctionSpaceTools.hpp"
00054 #include "Teuchos_oblackholestream.hpp"
00055 #include "Teuchos_RCP.hpp"
00056 #include "Teuchos_GlobalMPISession.hpp"
00057 
00058 using namespace std;
00059 using namespace Intrepid;
00060 
00061 #define INTREPID_TEST_COMMAND( S , throwCounter, nException )                                                              \
00062 {                                                                                                                          \
00063   ++nException;                                                                                                            \
00064   try {                                                                                                                    \
00065     S ;                                                                                                                    \
00066   }                                                                                                                        \
00067   catch (std::logic_error err) {                                                                                           \
00068       ++throwCounter;                                                                                                      \
00069       *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
00070       *outStream << err.what() << '\n';                                                                                    \
00071       *outStream << "-------------------------------------------------------------------------------" << "\n\n";           \
00072   };                                                                                                                       \
00073 }
00074 
00075 
00076 int main(int argc, char *argv[]) {
00077 
00078   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00079 
00080   // This little trick lets us print to std::cout only if
00081   // a (dummy) command-line argument is provided.
00082   int iprint     = argc - 1;
00083   Teuchos::RCP<std::ostream> outStream;
00084   Teuchos::oblackholestream bhs; // outputs nothing
00085   if (iprint > 0)
00086     outStream = Teuchos::rcp(&std::cout, false);
00087   else
00088     outStream = Teuchos::rcp(&bhs, false);
00089 
00090   // Save the format state of the original std::cout.
00091   Teuchos::oblackholestream oldFormatState;
00092   oldFormatState.copyfmt(std::cout);
00093 
00094   *outStream \
00095     << "===============================================================================\n" \
00096     << "|                                                                             |\n" \
00097     << "|               Unit Test (Basis_HGRAD_LINE_Cn_FEM_JACOBI)                    |\n" \
00098     << "|                                                                             |\n" \
00099     << "|     1) Conversion of Dof tags into Dof ordinals and back                    |\n" \
00100     << "|     2) Basis values for VALUE, GRAD, CURL, and Dk operators                 |\n" \
00101     << "|                                                                             |\n" \
00102     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00103     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00104     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00105     << "|                                                                             |\n" \
00106     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00107     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00108     << "|                                                                             |\n" \
00109     << "===============================================================================\n"\
00110     << "| TEST 1: Basis creation, exception testing                                   |\n"\
00111     << "===============================================================================\n";
00112 
00113   
00114   // Define basis and error flag
00115   double alpha = 0.0, beta = 0.0;
00116   Basis_HGRAD_LINE_Cn_FEM_JACOBI<double, FieldContainer<double> > lineBasis(5, alpha, beta);
00117   int errorFlag = 0;
00118 
00119   // Initialize throw counter for exception testing
00120   int nException     = 0;
00121   int throwCounter   = 0;
00122   
00123   // Define array containing vertices of the reference Line and a few other points   
00124   int numIntervals = 100;
00125   FieldContainer<double> lineNodes(numIntervals+1, 1);
00126   for (int i=0; i<numIntervals+1; i++) {
00127     lineNodes(i,0) = -1.0+(2.0*(double)i)/(double)numIntervals;
00128   }
00129 
00130   // Generic array for the output values; needs to be properly resized depending on the operator type
00131   FieldContainer<double> vals;
00132 
00133   try{
00134     // Exceptions 1-5: 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 #1
00137     INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(2,0,0), throwCounter, nException );
00138     // exception #2
00139     INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(1,1,1), throwCounter, nException );
00140     // exception #3
00141     INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(1,0,7), throwCounter, nException );
00142     // not an exception
00143     INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(1,0,5), throwCounter, nException ); --nException;
00144     // exception #4
00145     INTREPID_TEST_COMMAND( lineBasis.getDofTag(6), throwCounter, nException );
00146     // exception #5
00147     INTREPID_TEST_COMMAND( lineBasis.getDofTag(-1), throwCounter, nException );
00148     // not an exception
00149     INTREPID_TEST_COMMAND( lineBasis.getDofTag(5), throwCounter, nException ); --nException;
00150 #ifdef HAVE_INTREPID_DEBUG
00151     // Exceptions 6-16 test exception handling with incorrectly dimensioned input/output arrays
00152     // exception #6: input points array must be of rank-2
00153     FieldContainer<double> badPoints1(4, 5, 3);
00154     INTREPID_TEST_COMMAND( lineBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00155 
00156     // exception #7: dimension 1 in the input point array must equal space dimension of the cell
00157     FieldContainer<double> badPoints2(4, 3);
00158     INTREPID_TEST_COMMAND( lineBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00159 
00160     // exception #8: output values must be of rank-2 for OPERATOR_VALUE
00161     FieldContainer<double> badVals1(4, 3, 1);
00162     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals1, lineNodes, OPERATOR_VALUE), throwCounter, nException );
00163 
00164     // exception #9: output values must be of rank-3 for OPERATOR_GRAD
00165     FieldContainer<double> badVals2(4, 3);
00166     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_GRAD), throwCounter, nException );
00167 
00168     // exception #10: output values must be of rank-3 for OPERATOR_CURL
00169     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_CURL), throwCounter, nException );
00170 
00171     // exception #11: output values must be of rank-2 for OPERATOR_DIV
00172     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_DIV), throwCounter, nException );
00173 
00174     // exception #12: output values must be of rank-2 for OPERATOR_D1
00175     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_D1), throwCounter, nException );
00176 
00177     // exception #13: incorrect 0th dimension of output array (must equal number of basis functions)
00178     FieldContainer<double> badVals3(lineBasis.getCardinality() + 1, lineNodes.dimension(0));
00179     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals3, lineNodes, OPERATOR_VALUE), throwCounter, nException );
00180 
00181     // exception #14: incorrect 1st dimension of output array (must equal number of points)
00182     FieldContainer<double> badVals4(lineBasis.getCardinality(), lineNodes.dimension(0) + 1);
00183     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals4, lineNodes, OPERATOR_VALUE), throwCounter, nException );
00184 
00185     // exception #15: incorrect 2nd dimension of output array (must equal spatial dimension)
00186     FieldContainer<double> badVals5(lineBasis.getCardinality(), lineNodes.dimension(0), 2);
00187     INTREPID_TEST_COMMAND( lineBasis.getValues(badVals5, lineNodes, OPERATOR_GRAD), throwCounter, nException );
00188 
00189     // not an exception
00190     FieldContainer<double> goodVals2(lineBasis.getCardinality(), lineNodes.dimension(0));
00191     INTREPID_TEST_COMMAND( lineBasis.getValues(goodVals2, lineNodes, OPERATOR_VALUE), throwCounter, nException ); --nException;
00192 #endif
00193 
00194   }
00195   catch (std::logic_error err) {
00196     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00197     *outStream << err.what() << '\n';
00198     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00199     errorFlag = -1000;
00200   };
00201 
00202   // Check if number of thrown exceptions matches the one we expect
00203   if (throwCounter != nException) {
00204     errorFlag++;
00205     *outStream << std::setw(70) << "FAILURE! Incorrect number of exceptions." << "\n";
00206   }
00207 
00208 
00209   *outStream \
00210     << "\n"
00211     << "===============================================================================\n"\
00212     << "| TEST 3: orthogonality of basis functions                                    |\n"\
00213     << "===============================================================================\n";
00214 
00215   outStream -> precision(20);
00216 
00217   try {
00218 
00219     // Check orthogonality property for Legendre polynomials.
00220     int maxorder = 10;
00221 
00222     DefaultCubatureFactory<double>  cubFactory;                                   // create factory
00223     shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());   // create cell topology
00224     for (int ordi=0; ordi < maxorder; ordi++) {
00225       //create left basis
00226       Teuchos::RCP<Basis<double,FieldContainer<double> > > lineBasisLeft =
00227         Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM_JACOBI<double,FieldContainer<double> >(ordi) );
00228 
00229       for (int ordj=0; ordj < maxorder; ordj++) {
00230 
00231         //create right basis
00232         Teuchos::RCP<Basis<double,FieldContainer<double> > > lineBasisRight =
00233           Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM_JACOBI<double,FieldContainer<double> >(ordj) );
00234 
00235         // get cubature points and weights
00236         Teuchos::RCP<Cubature<double> > lineCub = cubFactory.create(line, ordi+ordj);
00237         int numPoints = lineCub->getNumPoints();
00238         FieldContainer<double> cubPoints (numPoints, lineCub->getDimension());
00239         FieldContainer<double> cubWeights(numPoints);
00240         FieldContainer<double> cubWeightsC(1, numPoints);
00241         lineCub->getCubature(cubPoints, cubWeights);
00242         // "reshape" weights
00243         for (int i=0; i<numPoints; i++) { cubWeightsC(0,i) = cubWeights(i); }
00244         
00245 
00246         // get basis values
00247         int numFieldsLeft  = lineBasisLeft ->getCardinality();
00248         int numFieldsRight = lineBasisRight->getCardinality();
00249         FieldContainer<double> valsLeft(numFieldsLeft,numPoints),
00250                                valsRight(numFieldsRight,numPoints);
00251         lineBasisLeft ->getValues(valsLeft,  cubPoints, OPERATOR_VALUE);
00252         lineBasisRight->getValues(valsRight, cubPoints, OPERATOR_VALUE);
00253 
00254         // reshape by cloning and integrate
00255         FieldContainer<double> valsLeftC(1, numFieldsLeft,numPoints),
00256                                valsRightC(1, numFieldsRight,numPoints),
00257                                massMatrix(1, numFieldsLeft, numFieldsRight);
00258         ArrayTools::cloneFields<double>(valsLeftC, valsLeft);
00259         ArrayTools::cloneFields<double>(valsRightC, valsRight);
00260         ArrayTools::scalarMultiplyDataField<double>(valsRightC, cubWeightsC, valsRightC);
00261         FunctionSpaceTools::integrate<double>(massMatrix, valsLeftC, valsRightC, COMP_CPP);
00262 
00263         // check orthogonality property
00264         for (int i=0; i<numFieldsLeft; i++) {
00265           for (int j=0; j<numFieldsRight; j++) {
00266 
00267             if (i==j) {
00268               if ( std::abs(massMatrix(0,i,j)-(double)(2.0/(2.0*j+1.0))) > INTREPID_TOL ) {
00269                 *outStream << "Incorrect ii (\"diagonal\") value for i=" << i << ", j=" << j << ": "
00270                            << massMatrix(0,i,j) << " != " << "2/(2*" << j << "+1)\n\n";
00271                 errorFlag++;
00272               }
00273             }
00274             else {
00275               if ( std::abs(massMatrix(0,i,j)) > INTREPID_TOL ) {
00276                 *outStream << "Incorrect ij (\"off-diagonal\") value for i=" << i << ", j=" << j << ": "
00277                            << massMatrix(0,i,j) << " != " << "0\n\n";
00278                 errorFlag++;
00279               }
00280             }
00281           }
00282         }
00283 
00284       }
00285     }
00286 
00287   }
00288   // Catch unexpected errors
00289   catch (std::logic_error err) {
00290     *outStream << err.what() << "\n\n";
00291     errorFlag = -1000;
00292   };
00293 
00294   *outStream \
00295     << "\n"
00296     << "===============================================================================\n"\
00297     << "| TEST 4: correctness of basis function derivatives                           |\n"\
00298     << "===============================================================================\n";
00299 
00300   outStream -> precision(20);
00301 
00302   // function values stored by bf, then pt
00303   double basisValues[] = {
00304     1.000000000000000, 1.000000000000000, 1.000000000000000,    \
00305     1.000000000000000, -1.000000000000000, -0.3333333333333333, \
00306     0.3333333333333333, 1.000000000000000, 1.000000000000000,   \
00307     -0.3333333333333333, -0.3333333333333333, 1.000000000000000,        \
00308     -1.000000000000000, 0.4074074074074074, -0.4074074074074074,        \
00309     1.000000000000000};
00310 
00311   double basisD1Values[] = 
00312     {0, 0, 0, 0, 1.000000000000000, 1.000000000000000, 1.000000000000000, \
00313      1.000000000000000, -3.000000000000000, -1.000000000000000,         \
00314      1.000000000000000, 3.000000000000000, 6.000000000000000,           \
00315      -0.6666666666666667, -0.6666666666666667, 6.000000000000000};
00316 
00317   double basisD2Values[] = 
00318     {0, 0, 0, 0, 0, 0, 0, 0, 3.000000000000000, 3.000000000000000,      \
00319      3.000000000000000, 3.000000000000000, -15.00000000000000,          \
00320      -5.000000000000000, 5.000000000000000, 15.00000000000000};
00321 
00322   double basisD3Values[] = 
00323     {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 15.00000000000000,             \
00324      15.00000000000000, 15.00000000000000, 15.00000000000000};
00325   
00326 
00327 
00328   try {
00329     Basis_HGRAD_LINE_Cn_FEM_JACOBI<double, FieldContainer<double> > lineBasis3(3, alpha, beta);
00330     int numIntervals = 3;
00331     FieldContainer<double> lineNodes3(numIntervals+1, 1);
00332     FieldContainer<double> vals;
00333     for (int i=0; i<numIntervals+1; i++) {
00334       lineNodes3(i,0) = -1.0+(2.0*(double)i)/(double)numIntervals;
00335     }
00336     int numFields = lineBasis3.getCardinality();
00337     int numPoints = lineNodes3.dimension(0);
00338 
00339     // test basis values
00340     vals.resize(numFields, numPoints);
00341     lineBasis3.getValues(vals,lineNodes3,OPERATOR_VALUE);
00342     for (int i = 0; i < numFields; i++) {
00343       for (int j = 0; j < numPoints; j++) {
00344         
00345         // Compute offset for (F,P) container
00346         int l =  j + i * numPoints;
00347         if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
00348           errorFlag++;
00349           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00350           
00351           // Output the multi-index of the value where the error is:
00352           *outStream << " At multi-index { ";
00353           *outStream << i << " ";*outStream << j << " ";
00354           *outStream << "}  computed value: " << vals(i,j)
00355                      << " but reference value: " << basisValues[l] << "\n";
00356         }
00357       }
00358     }
00359 
00360     // test basis derivatives
00361     vals.resize(numFields, numPoints,1);
00362     lineBasis3.getValues(vals,lineNodes3,OPERATOR_D1);
00363     for (int i = 0; i < numFields; i++) {
00364       for (int j = 0; j < numPoints; j++) {
00365         
00366         // Compute offset for (F,P) container
00367         int l =  j + i * numPoints;
00368         if (std::abs(vals(i,j,0) - basisD1Values[l]) > INTREPID_TOL) {
00369           errorFlag++;
00370           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00371           
00372           // Output the multi-index of the value where the error is:
00373           *outStream << " At multi-index { ";
00374           *outStream << i << " ";*outStream << j << " ";
00375           *outStream << "}  computed value: " << vals(i,j,0)
00376                      << " but reference value: " << basisD1Values[l] << "\n";
00377         }
00378       }
00379     }
00380 
00381     vals.resize(numFields, numPoints,1);
00382     lineBasis3.getValues(vals,lineNodes3,OPERATOR_D2);
00383     for (int i = 0; i < numFields; i++) {
00384       for (int j = 0; j < numPoints; j++) {
00385         
00386         // Compute offset for (F,P) container
00387         int l =  j + i * numPoints;
00388         if (std::abs(vals(i,j,0) - basisD2Values[l]) > INTREPID_TOL) {
00389           errorFlag++;
00390           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00391           
00392           // Output the multi-index of the value where the error is:
00393           *outStream << " At multi-index { ";
00394           *outStream << i << " ";*outStream << j << " ";
00395           *outStream << "}  computed value: " << vals(i,j,0)
00396                      << " but reference value: " << basisD2Values[l] << "\n";
00397         }
00398       }
00399     }
00400 
00401     vals.resize(numFields, numPoints,1);
00402     lineBasis3.getValues(vals,lineNodes3,OPERATOR_D3);
00403     for (int i = 0; i < numFields; i++) {
00404       for (int j = 0; j < numPoints; j++) {
00405         
00406         // Compute offset for (F,P) container
00407         int l =  j + i * numPoints;
00408         if (std::abs(vals(i,j,0) - basisD3Values[l]) > INTREPID_TOL) {
00409           errorFlag++;
00410           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00411           
00412           // Output the multi-index of the value where the error is:
00413           *outStream << " At multi-index { ";
00414           *outStream << i << " ";*outStream << j << " ";
00415           *outStream << "}  computed value: " << vals(i,j,0)
00416                      << " but reference value: " << basisD3Values[l] << "\n";
00417         }
00418       }
00419     }
00420   }
00421   // Catch unexpected errors
00422   catch (std::logic_error err) {
00423     *outStream << err.what() << "\n\n";
00424     errorFlag = -1000;
00425   };
00426 
00427 
00428   if (errorFlag != 0)
00429     std::cout << "End Result: TEST FAILED\n";
00430   else
00431     std::cout << "End Result: TEST PASSED\n";
00432 
00433   // reset format state of std::cout
00434   std::cout.copyfmt(oldFormatState);
00435 
00436   return errorFlag;
00437 }