Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Basis/HCURL_QUAD_In_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_HCURL_QUAD_In_FEM.hpp"
00037 #include "Intrepid_PointTools.hpp"
00038 #include "Teuchos_oblackholestream.hpp"
00039 #include "Teuchos_RCP.hpp"
00040 #include "Teuchos_GlobalMPISession.hpp"
00041 
00042 using namespace std;
00043 using namespace Intrepid;
00044 
00045 #define INTREPID_TEST_COMMAND( S , throwCounter, nException )                                                              \
00046 {                                                                                                                          \
00047   ++nException;                                                                                                            \
00048   try {                                                                                                                    \
00049     S ;                                                                                                                    \
00050   }                                                                                                                        \
00051   catch (std::logic_error err) {                                                                                           \
00052       ++throwCounter;                                                                                                      \
00053       *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
00054       *outStream << err.what() << '\n';                                                                                    \
00055       *outStream << "-------------------------------------------------------------------------------" << "\n\n";           \
00056   };                                                                                                                       \
00057 }
00058 
00059 int main(int argc, char *argv[]) {
00060 
00061   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00062 
00063   // This little trick lets us print to std::cout only if
00064   // a (dummy) command-line argument is provided.
00065   int iprint     = argc - 1;
00066   Teuchos::RCP<std::ostream> outStream;
00067   Teuchos::oblackholestream bhs; // outputs nothing
00068   if (iprint > 0)
00069     outStream = Teuchos::rcp(&std::cout, false);
00070   else
00071     outStream = Teuchos::rcp(&bhs, false);
00072   
00073   // Save the format state of the original std::cout.
00074   Teuchos::oblackholestream oldFormatState;
00075   oldFormatState.copyfmt(std::cout);
00076   
00077   *outStream \
00078     << "===============================================================================\n" \
00079     << "|                                                                             |\n" \
00080     << "|                 Unit Test (Basis_HCURL_QUAD_In_FEM)                          |\n" \
00081     << "|                                                                             |\n" \
00082     << "|     1) Conversion of Dof tags into Dof ordinals and back                    |\n" \
00083     << "|     2) Basis values for VALUE and CURL operators                            |\n" \
00084     << "|                                                                             |\n" \
00085     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00086     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00087     << "|                      Kara Peterson (kjpeter@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   // Define basis and error flag
00097   const int deg = 1;
00098   shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >()); 
00099 
00100   Basis_HCURL_QUAD_In_FEM<double, FieldContainer<double> > quadBasis(deg,POINTTYPE_EQUISPACED);
00101   int errorFlag = 0;
00102 
00103   // Initialize throw counter for exception testing
00104   int nException     = 0;
00105   int throwCounter   = 0;
00106 
00107   // Array with the 4 vertices of the reference Quadrilateral, its center and 4 more points
00108   FieldContainer<double> quadNodes(9, 2);
00109   quadNodes(0,0) = -1.0;  quadNodes(0,1) = -1.0;  
00110   quadNodes(1,0) =  0.0;  quadNodes(1,1) = -1.0;
00111   quadNodes(2,0) =  1.0;  quadNodes(2,1) = -1.0;
00112   quadNodes(3,0) = -1.0;  quadNodes(3,1) =  0.0;  
00113   quadNodes(4,0) =  0.0;  quadNodes(4,1) =  0.0;
00114   quadNodes(5,0) =  1.0;  quadNodes(5,1) =  0.0;
00115   quadNodes(6,0) = -1.0;  quadNodes(6,1) =  1.0;  
00116   quadNodes(7,0) =  0.0;  quadNodes(7,1) =  1.0;
00117   quadNodes(8,0) =  1.0;  quadNodes(8,1) =  1.0;
00118     
00119 
00120   // Generic array for the output values; needs to be properly resized depending on the operator type
00121   FieldContainer<double> vals;
00122 
00123   try{
00124     // exception #1: GRAD cannot be applied to HCURL functions 
00125     // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
00126     vals.resize(quadBasis.getCardinality(), quadNodes.dimension(0), 4 );
00127     INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_GRAD), throwCounter, nException );
00128 
00129     // exception #2: DIV cannot be applied to HCURL functions
00130     // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
00131     vals.resize(quadBasis.getCardinality(), quadNodes.dimension(0) );
00132     INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_DIV), throwCounter, nException );
00133         
00134     // Exceptions 3-7: 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 #3
00137     INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(3,0,0), throwCounter, nException );
00138     // exception #4
00139     INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(1,1,1), throwCounter, nException );
00140     // exception #5
00141     INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(0,4,1), throwCounter, nException );
00142     // exception #6
00143     INTREPID_TEST_COMMAND( quadBasis.getDofTag(12), throwCounter, nException );
00144     // exception #7
00145     INTREPID_TEST_COMMAND( quadBasis.getDofTag(-1), throwCounter, nException );
00146     
00147 #ifdef HAVE_INTREPID_DEBUG
00148     // Exceptions 8-15 test exception handling with incorrectly dimensioned input/output arrays
00149     // exception #8: input points array must be of rank-2
00150     FieldContainer<double> badPoints1(4, 5, 3);
00151     INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
00152     
00153     // exception #9 dimension 1 in the input point array must equal space dimension of the cell
00154     FieldContainer<double> badPoints2(4, quadBasis.getBaseCellTopology().getDimension() + 1);
00155     INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
00156     
00157     // exception #10 output values must be of rank-3 for OPERATOR_VALUE in 2D
00158     FieldContainer<double> badVals1(4, 3);
00159     INTREPID_TEST_COMMAND( quadBasis.getValues(badVals1, quadNodes, OPERATOR_VALUE), throwCounter, nException );
00160  
00161     FieldContainer<double> badCurls1(4,3,2);
00162     // exception #11 output values must be of rank-2 for OPERATOR_CURL
00163     INTREPID_TEST_COMMAND( quadBasis.getValues(badCurls1, quadNodes, OPERATOR_CURL), throwCounter, nException );
00164     
00165     // exception #12 incorrect 0th dimension of output array (must equal number of basis functions)
00166     FieldContainer<double> badVals2(quadBasis.getCardinality() + 1, quadNodes.dimension(0), quadBasis.getBaseCellTopology().getDimension());
00167     INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ;
00168     
00169     // exception #13 incorrect 1st  dimension of output array (must equal number of points)
00170     FieldContainer<double> badVals3(quadBasis.getCardinality(), quadNodes.dimension(0) + 1, quadBasis.getBaseCellTopology().getDimension() );
00171     INTREPID_TEST_COMMAND( quadBasis.getValues(badVals3, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ;
00172 
00173     // exception #14: incorrect 2nd dimension of output array for VALUE (must equal the space dimension)
00174     FieldContainer<double> badVals4(quadBasis.getCardinality(), quadNodes.dimension(0), quadBasis.getBaseCellTopology().getDimension() - 1);
00175     INTREPID_TEST_COMMAND( quadBasis.getValues(badVals4, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ;
00176         
00177     // exception #15: D2 cannot be applied to HCURL functions 
00178     // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
00179     vals.resize(quadBasis.getCardinality(), 
00180                 quadNodes.dimension(0),  
00181                 Intrepid::getDkCardinality(OPERATOR_D2, quadBasis.getBaseCellTopology().getDimension()));
00182     INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_D2), throwCounter, nException );
00183 #endif
00184     
00185   }
00186   catch (std::logic_error err) {
00187     *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
00188     *outStream << err.what() << '\n';
00189     *outStream << "-------------------------------------------------------------------------------" << "\n\n";
00190     errorFlag = -1000;
00191   };
00192   
00193   // Check if number of thrown exceptions matches the one we expect 
00194   // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match.
00195   if (throwCounter != nException) {
00196     errorFlag++;
00197     *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00198   }
00199 //#endif
00200   
00201   *outStream \
00202     << "\n"
00203     << "===============================================================================\n"\
00204     << "| TEST 2: correctness of tag to enum and enum to tag lookups                  |\n"\
00205     << "===============================================================================\n";
00206   
00207   try{
00208     std::vector<std::vector<int> > allTags = quadBasis.getAllDofTags();
00209     
00210     // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
00211     for (unsigned i = 0; i < allTags.size(); i++) {
00212       int bfOrd  = quadBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
00213       
00214       std::vector<int> myTag = quadBasis.getDofTag(bfOrd);
00215        if( !( (myTag[0] == allTags[i][0]) &&
00216               (myTag[1] == allTags[i][1]) &&
00217               (myTag[2] == allTags[i][2]) &&
00218               (myTag[3] == allTags[i][3]) ) ) {
00219         errorFlag++;
00220         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00221         *outStream << " getDofOrdinal( {" 
00222           << allTags[i][0] << ", " 
00223           << allTags[i][1] << ", " 
00224           << allTags[i][2] << ", " 
00225           << allTags[i][3] << "}) = " << bfOrd <<" but \n";   
00226         *outStream << " getDofTag(" << bfOrd << ") = { "
00227           << myTag[0] << ", " 
00228           << myTag[1] << ", " 
00229           << myTag[2] << ", " 
00230           << myTag[3] << "}\n";        
00231       }
00232     }
00233     
00234     // Now do the same but loop over basis functions
00235     for( int bfOrd = 0; bfOrd < quadBasis.getCardinality(); bfOrd++) {
00236       std::vector<int> myTag  = quadBasis.getDofTag(bfOrd);
00237       int myBfOrd = quadBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
00238       if( bfOrd != myBfOrd) {
00239         errorFlag++;
00240         *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00241         *outStream << " getDofTag(" << bfOrd << ") = { "
00242           << myTag[0] << ", " 
00243           << myTag[1] << ", " 
00244           << myTag[2] << ", " 
00245           << myTag[3] << "} but getDofOrdinal({" 
00246           << myTag[0] << ", " 
00247           << myTag[1] << ", " 
00248           << myTag[2] << ", " 
00249           << myTag[3] << "} ) = " << myBfOrd << "\n";
00250       }
00251     }
00252   }
00253   catch (std::logic_error err){
00254     *outStream << err.what() << "\n\n";
00255     errorFlag = -1000;
00256   };
00257   
00258   *outStream \
00259     << "\n"
00260     << "===============================================================================\n"\
00261     << "| TEST 3: correctness of basis function values                                |\n"\
00262     << "===============================================================================\n";
00263   
00264   outStream -> precision(20);
00265   
00266   // VALUE: Each row pair gives the 4x2 correct basis set values at an evaluation point: (P,F,D) layout
00267   double basisValues[] = {
00268     // first bf, first row of points
00269     0.0, 1.0, 0.0, 0.5, 0.0, 0.0, 
00270     // first bf, second row of points
00271     0.0, 1.0, 0.0, 0.5, 0.0, 0.0, 
00272     // first bf, third row of points
00273     0.0, 1.0, 0.0, 0.5, 0.0, 0.0, 
00274     // second bf, first row of points
00275     0.0, 0.0, 0.0, 0.5, 0.0, 1.0, 
00276     // second bf, second row of points
00277     0.0, 0.0, 0.0, 0.5, 0.0, 1.0, 
00278     // second bf, third row of points
00279     0.0, 0.0, 0.0, 0.5, 0.0, 1.0, 
00280     // third bf, first row of points
00281     1.0, 0.0, 1.0, 0.0, 1.0, 0.0,
00282     // third bf, second row of points
00283     0.5, 0.0, 0.5, 0.0, 0.5, 0.0,
00284     // third bf, third row of points
00285     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 
00286     // fourth bf, first row of points
00287     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 
00288     // fourth bf, second row of points
00289     0.5, 0.0, 0.5, 0.0, 0.5, 0.0,
00290     // fourth bf, third row of points
00291     1.0, 0.0, 1.0, 0.0, 1.0, 0.0
00292   };
00293   
00294   // CURL: correct values in (F,P) format
00295   double basisCurls[] = {
00296     -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
00297     0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
00298     -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
00299     0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5
00300   };
00301 
00302   
00303   try{
00304         
00305     // Dimensions for the output arrays:
00306     int numFields = quadBasis.getCardinality();
00307     int numPoints = quadNodes.dimension(0);
00308     int spaceDim  = quadBasis.getBaseCellTopology().getDimension();
00309     
00310     // Generic array for values and curls that will be properly sized before each call
00311     FieldContainer<double> vals;
00312     
00313     // Check VALUE of basis functions: resize vals to rank-3 container:
00314     vals.resize(numFields, numPoints, spaceDim);
00315     quadBasis.getValues(vals, quadNodes, OPERATOR_VALUE);
00316     for (int i = 0; i < numFields; i++) {
00317       for (int j = 0; j < numPoints; j++) {
00318         for (int k = 0; k < spaceDim; k++) {
00319            
00320           // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k
00321           int l = k + i * spaceDim * numPoints + j * spaceDim;
00322            if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
00323              errorFlag++;
00324              *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00325         
00326              // Output the multi-index of the value where the error is: 
00327              *outStream << " At multi-index { ";
00328              *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
00329              *outStream << "}  computed value: " << vals(i,j,k)
00330                << " but reference value: " << basisValues[l] << "\n";
00331             }
00332          }
00333       }
00334     }
00335 
00336     // Check CURL of basis function: resize vals to rank-2 container
00337     vals.resize(numFields, numPoints);
00338     quadBasis.getValues(vals, quadNodes, OPERATOR_CURL);
00339     for (int i = 0; i < numFields; i++) {
00340       for (int j = 0; j < numPoints; j++) {
00341         int l =  j + i * numPoints;
00342         if (std::abs(vals(i,j) - basisCurls[l]) > INTREPID_TOL) {
00343           errorFlag++;
00344           *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
00345           
00346           // Output the multi-index of the value where the error is:
00347           *outStream << " At multi-index { ";
00348           *outStream << i << " ";*outStream << j << " ";
00349           *outStream << "}  computed curl component: " << vals(i,j)
00350             << " but reference curl component: " << basisCurls[l] << "\n";
00351         }
00352       }
00353     }
00354     
00355    }    
00356   
00357   // Catch unexpected errors
00358   catch (std::logic_error err) {
00359     *outStream << err.what() << "\n\n";
00360     errorFlag = -1000;
00361   };
00362   
00363   if (errorFlag != 0)
00364     std::cout << "End Result: TEST FAILED\n";
00365   else
00366     std::cout << "End Result: TEST PASSED\n";
00367   
00368   // reset format state of std::cout
00369   std::cout.copyfmt(oldFormatState);
00370   
00371   return errorFlag;
00372 }