Intrepid
http://trilinos.sandia.gov/packages/docs/r10.12/packages/intrepid/example/Drivers/example_03NL.cpp
Go to the documentation of this file.
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 
00074 // Intrepid includes
00075 #include "Intrepid_FunctionSpaceTools.hpp"
00076 #include "Intrepid_FieldContainer.hpp"
00077 #include "Intrepid_CellTools.hpp"
00078 #include "Intrepid_ArrayTools.hpp"
00079 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
00080 #include "Intrepid_RealSpaceTools.hpp"
00081 #include "Intrepid_DefaultCubatureFactory.hpp"
00082 #include "Intrepid_Utils.hpp"
00083 
00084 // Epetra includes
00085 #include "Epetra_Time.h"
00086 #include "Epetra_Map.h"
00087 #include "Epetra_FECrsMatrix.h"
00088 #include "Epetra_FEVector.h"
00089 #include "Epetra_SerialComm.h"
00090 
00091 // Teuchos includes
00092 #include "Teuchos_oblackholestream.hpp"
00093 #include "Teuchos_RCP.hpp"
00094 #include "Teuchos_BLAS.hpp"
00095 #include "Teuchos_Time.hpp"
00096 
00097 // Shards includes
00098 #include "Shards_CellTopology.hpp"
00099 
00100 // EpetraExt includes
00101 #include "EpetraExt_RowMatrixOut.h"
00102 #include "EpetraExt_MultiVectorOut.h"
00103 #include "EpetraExt_MatrixMatrix.h"
00104 
00105 // Sacado includes
00106 #include "Sacado.hpp"
00107 #include "Sacado_Fad_DVFad.hpp"
00108 #include "Sacado_Fad_SimpleFad.hpp"
00109 #include "Sacado_CacheFad_DFad.hpp"
00110 #include "Sacado_CacheFad_SFad.hpp"
00111 #include "Sacado_CacheFad_SLFad.hpp"
00112 
00113 
00114 using namespace std;
00115 using namespace Intrepid;
00116 
00117 #define INTREPID_INTEGRATE_COMP_ENGINE COMP_BLAS
00118 
00119 #define BATCH_SIZE 10
00120 
00121 //typedef Sacado::Fad::DFad<double> FadType;
00122 //typedef Sacado::CacheFad::DFad<double> FadType;
00123 //typedef Sacado::ELRCacheFad::DFad<double> FadType;
00124 //typedef Sacado::Fad::SFad<double,8> FadType;
00125 typedef Sacado::CacheFad::SFad<double,8> FadType;
00126 //typedef Sacado::ELRCacheFad::SFad<double,8> FadType;
00127 //typedef Sacado::Fad::SLFad<double,8> FadType;
00128 //typedef Sacado::CacheFad::SLFad<double,8> FadType;
00129 //typedef Sacado::ELRCacheFad::SLFad<double,8> FadType;
00130 
00131 //#define DUMP_DATA
00132 
00133 // Functions to evaluate nonlinear terms
00134 void dfunc_u(FieldContainer<double>, FieldContainer<double>);
00135 
00136 template<class ScalarT>
00137 void func_u(FieldContainer<ScalarT>, FieldContainer<ScalarT>);
00138 // 
00139 
00140 int main(int argc, char *argv[]) {
00141 
00142     // Check number of arguments
00143     if (argc < 4) {
00144       std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n";
00145       std::cout <<"Usage:\n\n";
00146       std::cout <<"  ./Intrepid_example_Drivers_Example_03NL.exe NX NY NZ verbose\n\n";
00147       std::cout <<" where \n";
00148       std::cout <<"   int NX              - num intervals in x direction (assumed box domain, 0,1) \n";
00149       std::cout <<"   int NY              - num intervals in y direction (assumed box domain, 0,1) \n";
00150       std::cout <<"   int NZ              - num intervals in z direction (assumed box domain, 0,1) \n";
00151       std::cout <<"   verbose (optional)  - any character, indicates verbose output \n\n";
00152       exit(1);
00153     }
00154   
00155     // This little trick lets us print to std::cout only if
00156     // a (dummy) command-line argument is provided.
00157     int iprint     = argc - 1;
00158     Teuchos::RCP<std::ostream> outStream;
00159     Teuchos::oblackholestream bhs; // outputs nothing
00160     if (iprint > 3)
00161       outStream = Teuchos::rcp(&std::cout, false);
00162     else
00163       outStream = Teuchos::rcp(&bhs, false);
00164   
00165     // Save the format state of the original std::cout.
00166     Teuchos::oblackholestream oldFormatState;
00167     oldFormatState.copyfmt(std::cout);
00168   
00169     *outStream \
00170     << "===============================================================================\n" \
00171     << "|                                                                             |\n" \
00172     << "|  Example: Generate PDE Jacobian for a Nonlinear Reaction-Diffusion          |\n" \
00173     << "|                   Equation on Hexahedral Mesh                               |\n" \
00174     << "|                                                                             |\n" \
00175     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00176     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00177     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00178     << "|                                                                             |\n" \
00179     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00180     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00181     << "|                                                                             |\n" \
00182     << "===============================================================================\n";
00183 
00184 
00185     // ************************************ GET INPUTS **************************************
00186 
00187     int NX = atoi(argv[1]);  // num intervals in x direction (assumed box domain, 0,1)
00188     int NY = atoi(argv[2]);  // num intervals in y direction (assumed box domain, 0,1)
00189     int NZ = atoi(argv[3]);  // num intervals in z direction (assumed box domain, 0,1)
00190 
00191     // *********************************** CELL TOPOLOGY **********************************
00192 
00193     // Get cell topology for base hexahedron
00194     typedef shards::CellTopology CellTopology;
00195     CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
00196 
00197     // Get dimensions 
00198     int numNodesPerElem = hex_8.getNodeCount();
00199     int spaceDim = hex_8.getDimension();
00200 
00201     // *********************************** GENERATE MESH ************************************
00202 
00203     *outStream << "Generating mesh ... \n\n";
00204 
00205     *outStream << "   NX" << "   NY" << "   NZ\n";
00206     *outStream << std::setw(5) << NX <<
00207                   std::setw(5) << NY <<
00208                   std::setw(5) << NZ << "\n\n";
00209 
00210     // Print mesh information
00211     int numElems = NX*NY*NZ;
00212     int numNodes = (NX+1)*(NY+1)*(NZ+1);
00213     *outStream << " Number of Elements: " << numElems << " \n";
00214     *outStream << "    Number of Nodes: " << numNodes << " \n\n";
00215 
00216     // Cube
00217     double leftX = 0.0, rightX = 1.0;
00218     double leftY = 0.0, rightY = 1.0;
00219     double leftZ = 0.0, rightZ = 1.0;
00220 
00221     // Mesh spacing
00222     double hx = (rightX-leftX)/((double)NX);
00223     double hy = (rightY-leftY)/((double)NY);
00224     double hz = (rightZ-leftZ)/((double)NZ);
00225 
00226     // Get nodal coordinates
00227     FieldContainer<double> nodeCoord(numNodes, spaceDim);
00228     FieldContainer<int> nodeOnBoundary(numNodes);
00229     int inode = 0;
00230     for (int k=0; k<NZ+1; k++) {
00231       for (int j=0; j<NY+1; j++) {
00232         for (int i=0; i<NX+1; i++) {
00233           nodeCoord(inode,0) = leftX + (double)i*hx;
00234           nodeCoord(inode,1) = leftY + (double)j*hy;
00235           nodeCoord(inode,2) = leftZ + (double)k*hz;
00236           if (k==0 || j==0 || i==0 || k==NZ || j==NY || i==NX){
00237              nodeOnBoundary(inode)=1;
00238           }
00239           else {
00240              nodeOnBoundary(inode)=0;
00241           }
00242           inode++;
00243         }
00244       }
00245     }
00246 
00247 #ifdef DUMP_DATA
00248     // Print nodal coords
00249     ofstream fcoordout("coords.dat");
00250     for (int i=0; i<numNodes; i++) {
00251        fcoordout << nodeCoord(i,0) <<" ";
00252        fcoordout << nodeCoord(i,1) <<" ";
00253        fcoordout << nodeCoord(i,2) <<"\n";
00254     }
00255     fcoordout.close();
00256 #endif
00257 
00258 
00259     // Element to Node map
00260     FieldContainer<int> elemToNode(numElems, numNodesPerElem);
00261     int ielem = 0;
00262     for (int k=0; k<NZ; k++) {
00263       for (int j=0; j<NY; j++) {
00264         for (int i=0; i<NX; i++) {
00265           elemToNode(ielem,0) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i;
00266           elemToNode(ielem,1) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i + 1;
00267           elemToNode(ielem,2) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i + 1;
00268           elemToNode(ielem,3) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i;
00269           elemToNode(ielem,4) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i;
00270           elemToNode(ielem,5) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i + 1;
00271           elemToNode(ielem,6) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i + 1;
00272           elemToNode(ielem,7) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i;
00273           ielem++;
00274         }
00275       }
00276     }
00277 #ifdef DUMP_DATA
00278     // Output connectivity
00279     ofstream fe2nout("elem2node.dat");
00280     for (int k=0; k<NZ; k++) {
00281       for (int j=0; j<NY; j++) {
00282         for (int i=0; i<NX; i++) {
00283           int ielem = i + j * NX + k * NX * NY;
00284           for (int m=0; m<numNodesPerElem; m++){
00285               fe2nout << elemToNode(ielem,m) <<"  ";
00286            }
00287           fe2nout <<"\n";
00288         }
00289       }
00290     }
00291     fe2nout.close();
00292 #endif
00293 
00294 
00295     // ************************************ CUBATURE **************************************
00296 
00297     *outStream << "Getting cubature ... \n\n";
00298 
00299     // Get numerical integration points and weights
00300     DefaultCubatureFactory<double>  cubFactory;                                   
00301     int cubDegree = 2;
00302     Teuchos::RCP<Cubature<double> > hexCub = cubFactory.create(hex_8, cubDegree); 
00303 
00304     int cubDim       = hexCub->getDimension();
00305     int numCubPoints = hexCub->getNumPoints();
00306 
00307     FieldContainer<double> cubPoints(numCubPoints, cubDim);
00308     FieldContainer<double> cubWeights(numCubPoints);
00309 
00310     hexCub->getCubature(cubPoints, cubWeights);
00311 
00312 
00313     // ************************************** BASIS ***************************************
00314 
00315     *outStream << "Getting basis ... \n\n";
00316 
00317     // Define basis 
00318     Basis_HGRAD_HEX_C1_FEM<double, FieldContainer<double> > hexHGradBasis;
00319     int numFieldsG = hexHGradBasis.getCardinality();
00320     FieldContainer<double> hexGVals(numFieldsG, numCubPoints); 
00321     FieldContainer<double> hexGrads(numFieldsG, numCubPoints, spaceDim); 
00322 
00323     // Evaluate basis values and gradients at cubature points
00324     hexHGradBasis.getValues(hexGVals, cubPoints, OPERATOR_VALUE);
00325     hexHGradBasis.getValues(hexGrads, cubPoints, OPERATOR_GRAD);
00326 
00327 
00328     // ******** FEM ASSEMBLY *************
00329 
00330     *outStream << "Building PDE Jacobian ... \n\n";
00331 
00332     // Settings and data structures for mass and stiffness matrices
00333     typedef CellTools<double>  CellTools;
00334     typedef FunctionSpaceTools fst;
00335     int numCells = BATCH_SIZE; 
00336     int numBatches = numElems/numCells; 
00337 
00338     // Container for nodes
00339     FieldContainer<double> hexNodes(numCells, numNodesPerElem, spaceDim);
00340     // Containers for Jacobian
00341     FieldContainer<double> hexJacobian(numCells, numCubPoints, spaceDim, spaceDim);
00342     FieldContainer<double> hexJacobInv(numCells, numCubPoints, spaceDim, spaceDim);
00343     FieldContainer<double> hexJacobDet(numCells, numCubPoints);
00344     // Containers for HGRAD bases
00345     FieldContainer<double> localPDEjacobian(numCells, numFieldsG, numFieldsG);
00346     FieldContainer<double> weightedMeasure(numCells, numCubPoints);
00347     FieldContainer<double> hexGValsTransformed(numCells, numFieldsG, numCubPoints);
00348     FieldContainer<double> hexGValsTransformedWeighted(numCells, numFieldsG, numCubPoints);
00349     FieldContainer<double> hexGradsTransformed(numCells, numFieldsG, numCubPoints, spaceDim);
00350     FieldContainer<double> hexGradsTransformedWeighted(numCells, numFieldsG, numCubPoints, spaceDim);
00351 
00352     // Global arrays in Epetra format 
00353     Epetra_SerialComm Comm;
00354     Epetra_Map globalMapG(numNodes, 0, Comm);
00355     Epetra_FECrsMatrix StiffMatrix(Copy, globalMapG, 64);
00356 
00357     // Additional arrays used in analytic assembly
00358     FieldContainer<double> u_coeffs(numCells, numFieldsG);
00359     FieldContainer<double> u_FE_val(numCells, numCubPoints);
00360     FieldContainer<double> df_of_u(numCells, numCubPoints);
00361     FieldContainer<double> df_of_u_times_basis(numCells, numFieldsG, numCubPoints);
00362 
00363 
00364     // Additional arrays used in AD-based assembly.
00365     FieldContainer<FadType> u_coeffsAD(numCells, numFieldsG);  
00366     FieldContainer<FadType> u_FE_gradAD(numCells, numCubPoints, spaceDim);
00367     FieldContainer<FadType> u_FE_valAD(numCells, numCubPoints);
00368     FieldContainer<FadType> f_of_u_AD(numCells, numCubPoints);
00369     FieldContainer<FadType> cellResidualAD(numCells, numFieldsG);
00370     for (int c=0; c<numCells; c++) {
00371       for(int f=0; f<numFieldsG; f++) {
00372           u_coeffsAD(c,f) = FadType(numFieldsG, f, 1.3);
00373       }
00374     }
00375 
00376     Teuchos::Time timer_jac_analytic("Time to compute element PDE Jacobians analytically: ");
00377     Teuchos::Time timer_jac_fad     ("Time to compute element PDE Jacobians using AD:     ");
00378     Teuchos::Time timer_jac_insert  ("Time for global insert,  w/o graph: ");
00379     Teuchos::Time timer_jac_insert_g("Time for global insert,  w/  graph: ");
00380     Teuchos::Time timer_jac_ga      ("Time for GlobalAssemble, w/o graph: ");
00381     Teuchos::Time timer_jac_ga_g    ("Time for GlobalAssemble, w/  graph: ");
00382     Teuchos::Time timer_jac_fc      ("Time for FillComplete,   w/o graph: ");
00383     Teuchos::Time timer_jac_fc_g    ("Time for FillComplete,   w/  graph: ");
00384 
00385 
00386 
00387 
00388     // *** Analytic element loop ***
00389     for (int bi=0; bi<numBatches; bi++) {
00390 
00391       // Physical cell coordinates
00392       for (int ci=0; ci<numCells; ci++) {
00393         int k = bi*numCells+ci;
00394         for (int i=0; i<numNodesPerElem; i++) {
00395             hexNodes(ci,i,0) = nodeCoord(elemToNode(k,i),0);
00396             hexNodes(ci,i,1) = nodeCoord(elemToNode(k,i),1);
00397             hexNodes(ci,i,2) = nodeCoord(elemToNode(k,i),2);
00398         }
00399       }
00400 
00401       // Compute cell Jacobians, their inverses and their determinants
00402       CellTools::setJacobian(hexJacobian, cubPoints, hexNodes, hex_8);
00403       CellTools::setJacobianInv(hexJacobInv, hexJacobian );
00404       CellTools::setJacobianDet(hexJacobDet, hexJacobian );
00405 
00406       // ******************** COMPUTE ELEMENT HGrad STIFFNESS MATRICES WITHOUT AD *******************
00407 
00408       // transform to physical coordinates 
00409       fst::HGRADtransformGRAD<double>(hexGradsTransformed, hexJacobInv, hexGrads);
00410       
00411       // compute weighted measure
00412       fst::computeCellMeasure<double>(weightedMeasure, hexJacobDet, cubWeights);
00413 
00414       // multiply values with weighted measure
00415       fst::multiplyMeasure<double>(hexGradsTransformedWeighted,
00416                                    weightedMeasure, hexGradsTransformed);
00417 
00418       // u_coeffs equals the value of u_coeffsAD
00419       for(int i=0; i<numFieldsG; i++){
00420         u_coeffs(0,i) = u_coeffsAD(0,i).val();
00421       }
00422 
00423       timer_jac_analytic.start(); // START TIMER
00424       // integrate to account for linear stiffness term
00425       fst::integrate<double>(localPDEjacobian, hexGradsTransformed, hexGradsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE);
00426 
00427       // represent value of the current state (iterate) as a linear combination of the basis functions
00428       u_FE_val.initialize();
00429       fst::evaluate<double>(u_FE_val, u_coeffs, hexGValsTransformed);
00430      
00431       // evaluate derivative of the nonlinear term and multiply by basis function
00432       dfunc_u(df_of_u, u_FE_val);
00433       fst::scalarMultiplyDataField<double>(df_of_u_times_basis, df_of_u, hexGValsTransformed);
00434 
00435       // integrate to account for nonlinear reaction term
00436       fst::integrate<double>(localPDEjacobian, df_of_u_times_basis, hexGValsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE, true);
00437       timer_jac_analytic.stop(); // STOP TIMER
00438 
00439       // assemble into global matrix
00440       for (int ci=0; ci<numCells; ci++) {
00441         int k = bi*numCells+ci;
00442         std::vector<int> rowIndex(numFieldsG);
00443         std::vector<int> colIndex(numFieldsG);
00444         for (int row = 0; row < numFieldsG; row++){
00445           rowIndex[row] = elemToNode(k,row);
00446         }
00447         for (int col = 0; col < numFieldsG; col++){
00448           colIndex[col] = elemToNode(k,col);
00449         }
00450         // We can insert an entire matrix at a time, but we opt for rows only.
00451         //timer_jac_insert.start();
00452         //StiffMatrix.InsertGlobalValues(numFieldsG, &rowIndex[0], numFieldsG, &colIndex[0], &localPDEjacobian(ci,0,0));
00453         //timer_jac_insert.stop();
00454         for (int row = 0; row < numFieldsG; row++){
00455           timer_jac_insert.start();
00456           StiffMatrix.InsertGlobalValues(1, &rowIndex[row], numFieldsG, &colIndex[0], &localPDEjacobian(ci,row,0));
00457           timer_jac_insert.stop();
00458         }
00459       }
00460 
00461     } // *** end analytic element loop ***
00462      
00463     // Assemble global objects
00464     timer_jac_ga.start(); StiffMatrix.GlobalAssemble(); timer_jac_ga.stop();
00465     timer_jac_fc.start(); StiffMatrix.FillComplete(); timer_jac_fc.stop();
00466 
00467 
00468 
00469 
00470     // *** AD element loop ***
00471 
00472     Epetra_CrsGraph mgraph = StiffMatrix.Graph();
00473     Epetra_FECrsMatrix StiffMatrixViaAD(Copy, mgraph);
00474 
00475     for (int bi=0; bi<numBatches; bi++) {
00476 
00477       // ******************** COMPUTE ELEMENT HGrad STIFFNESS MATRICES AND RIGHT-HAND SIDE WITH AD ********************
00478 
00479       // Physical cell coordinates
00480       for (int ci=0; ci<numCells; ci++) {
00481         int k = bi*numCells+ci;
00482         for (int i=0; i<numNodesPerElem; i++) {
00483             hexNodes(ci,i,0) = nodeCoord(elemToNode(k,i),0);
00484             hexNodes(ci,i,1) = nodeCoord(elemToNode(k,i),1);
00485             hexNodes(ci,i,2) = nodeCoord(elemToNode(k,i),2);
00486         }
00487       }
00488 
00489       // Compute cell Jacobians, their inverses and their determinants
00490       CellTools::setJacobian(hexJacobian, cubPoints, hexNodes, hex_8);
00491       CellTools::setJacobianInv(hexJacobInv, hexJacobian );
00492       CellTools::setJacobianDet(hexJacobDet, hexJacobian );
00493 
00494       // transform to physical coordinates
00495       fst::HGRADtransformGRAD<double>(hexGradsTransformed, hexJacobInv, hexGrads);
00496     
00497       // compute weighted measure
00498       fst::computeCellMeasure<double>(weightedMeasure, hexJacobDet, cubWeights);
00499     
00500       // multiply values with weighted measure
00501       fst::multiplyMeasure<double>(hexGradsTransformedWeighted, weightedMeasure, hexGradsTransformed);
00502 
00503       // transform basis values to physical coordinates 
00504       fst::HGRADtransformVALUE<double>(hexGValsTransformed, hexGVals);
00505 
00506       // multiply values with weighted measure
00507       fst::multiplyMeasure<double>(hexGValsTransformedWeighted,
00508                                    weightedMeasure, hexGValsTransformed);
00509 
00510       timer_jac_fad.start(); // START TIMER
00511       // represent gradient of the current state (iterate) as a linear combination of the gradients of basis functions
00512       // use AD arrays !
00513       u_FE_gradAD.initialize();
00514       fst::evaluate<FadType>(u_FE_gradAD, u_coeffsAD, hexGradsTransformed);
00515 
00516       // represent value of the current state (iterate) as a linear combination of the basis functions
00517       // use AD arrays !
00518       u_FE_valAD.initialize();
00519       fst::evaluate<FadType>(u_FE_valAD, u_coeffsAD, hexGValsTransformed);
00520       // compute nonlinear term
00521       func_u(f_of_u_AD, u_FE_valAD);
00522 
00523       // integrate to compute element residual   
00524       fst::integrate<FadType>(cellResidualAD, u_FE_gradAD,  hexGradsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE);
00525       fst::integrate<FadType>(cellResidualAD, f_of_u_AD, hexGValsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE, true);
00526       timer_jac_fad.stop(); // STOP TIMER
00527 
00528       // assemble into global matrix
00529       for (int ci=0; ci<numCells; ci++) {
00530         int k = bi*numCells+ci;
00531         std::vector<int> rowIndex(numFieldsG);
00532         std::vector<int> colIndex(numFieldsG);
00533         for (int row = 0; row < numFieldsG; row++){
00534           rowIndex[row] = elemToNode(k,row);
00535         }
00536         for (int col = 0; col < numFieldsG; col++){
00537           colIndex[col] = elemToNode(k,col);
00538         }
00539         for (int row = 0; row < numFieldsG; row++){
00540           timer_jac_insert_g.start();
00541           StiffMatrixViaAD.SumIntoGlobalValues(1, &rowIndex[row], numFieldsG, &colIndex[0], cellResidualAD(ci,row).dx());
00542           timer_jac_insert_g.stop();
00543         }
00544       }
00545  
00546     } // *** end AD element loop ***
00547 
00548     // Assemble global objects
00549     timer_jac_ga_g.start(); StiffMatrixViaAD.GlobalAssemble(); timer_jac_ga_g.stop();
00550     timer_jac_fc_g.start(); StiffMatrixViaAD.FillComplete(); timer_jac_fc_g.stop();
00551 
00552 
00553 
00554     /****** Output *******/
00555 
00556 #ifdef DUMP_DATA
00557     // Dump matrices to disk
00558     EpetraExt::RowMatrixToMatlabFile("stiff_matrix.dat",StiffMatrix);
00559     EpetraExt::RowMatrixToMatlabFile("stiff_matrixAD.dat",StiffMatrixViaAD);
00560 #endif
00561 
00562     // take the infinity norm of the difference between StiffMatrix and StiffMatrixViaAD to see that 
00563     // the two matrices are the same
00564     EpetraExt::MatrixMatrix::Add(StiffMatrix, false, 1.0, StiffMatrixViaAD, -1.0);
00565     double normMat = StiffMatrixViaAD.NormInf();
00566     *outStream << "Infinity norm of difference between stiffness matrices = " << normMat << "\n";
00567 
00568 
00569     *outStream << "\n\nNumber of global nonzeros: " << StiffMatrix.NumGlobalNonzeros() << "\n\n";
00570 
00571     *outStream << timer_jac_analytic.name() << " " << timer_jac_analytic.totalElapsedTime() << " sec\n";
00572     *outStream << timer_jac_fad.name()      << " " << timer_jac_fad.totalElapsedTime()      << " sec\n\n";
00573     *outStream << timer_jac_insert.name()   << " " << timer_jac_insert.totalElapsedTime()   << " sec\n";
00574     *outStream << timer_jac_insert_g.name() << " " << timer_jac_insert_g.totalElapsedTime() << " sec\n\n";
00575     *outStream << timer_jac_ga.name()       << " " << timer_jac_ga.totalElapsedTime()       << " sec\n";
00576     *outStream << timer_jac_ga_g.name()     << " " << timer_jac_ga_g.totalElapsedTime()     << " sec\n\n";
00577     *outStream << timer_jac_fc.name()       << " " << timer_jac_fc.totalElapsedTime()       << " sec\n";
00578     *outStream << timer_jac_fc_g.name()     << " " << timer_jac_fc_g.totalElapsedTime()     << " sec\n\n";
00579 
00580     if ((normMat < 1.0e4*INTREPID_TOL)) {
00581       std::cout << "End Result: TEST PASSED\n";
00582     }
00583     else {
00584       std::cout << "End Result: TEST FAILED\n";
00585     }
00586    
00587     // reset format state of std::cout
00588     std::cout.copyfmt(oldFormatState);
00589    
00590     return 0;
00591 }
00592 
00593 
00594 template<class ScalarT>
00595 void func_u(FieldContainer<ScalarT> fu, FieldContainer<ScalarT> u) {
00596   int num_cells = u.dimension(0);
00597   int num_cub_p = u.dimension(1);
00598   for(int c=0; c<num_cells; c++){
00599     for(int p=0; p<num_cub_p; p++){
00600       fu(c,p) = std::pow(u(c,p),3) + std::exp(u(c,p));
00601     }
00602   }
00603 }
00604 
00605 
00606 void dfunc_u(FieldContainer<double> dfu, FieldContainer<double> u) {
00607   int num_cells = u.dimension(0);
00608   int num_cub_p = u.dimension(1);
00609   for(int c=0; c<num_cells; c++) {
00610     for(int p=0; p<num_cub_p; p++) {
00611       dfu(c,p) = 3*u(c,p)*u(c,p) + std::exp(u(c,p));
00612     }
00613   }
00614 }