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