Intrepid
example_03AD.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Pavel Bochev (pbboche@sandia.gov)
38 // Denis Ridzal (dridzal@sandia.gov), or
39 // Kara Peterson (kjpeter@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
81 // Intrepid includes
84 #include "Intrepid_CellTools.hpp"
85 #include "Intrepid_ArrayTools.hpp"
86 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
89 #include "Intrepid_Utils.hpp"
90 
91 // Epetra includes
92 #include "Epetra_Time.h"
93 #include "Epetra_Map.h"
94 #include "Epetra_FECrsMatrix.h"
95 #include "Epetra_FEVector.h"
96 #include "Epetra_SerialComm.h"
97 
98 // Teuchos includes
99 #include "Teuchos_oblackholestream.hpp"
100 #include "Teuchos_RCP.hpp"
101 #include "Teuchos_BLAS.hpp"
102 #include "Teuchos_Time.hpp"
103 
104 // Shards includes
105 #include "Shards_CellTopology.hpp"
106 
107 // EpetraExt includes
108 #include "EpetraExt_RowMatrixOut.h"
109 #include "EpetraExt_MultiVectorOut.h"
110 #include "EpetraExt_MatrixMatrix.h"
111 
112 // Sacado includes
113 #include "Sacado.hpp"
114 #include "Sacado_Fad_DVFad.hpp"
115 #include "Sacado_Fad_SimpleFad.hpp"
116 #include "Sacado_CacheFad_DFad.hpp"
117 #include "Sacado_CacheFad_SFad.hpp"
118 #include "Sacado_CacheFad_SLFad.hpp"
119 
120 
121 using namespace std;
122 using namespace Intrepid;
123 
124 #define INTREPID_INTEGRATE_COMP_ENGINE COMP_BLAS
125 
126 #define BATCH_SIZE 10
127 
128 //typedef Sacado::Fad::DFad<double> FadType;
129 //typedef Sacado::CacheFad::DFad<double> FadType;
130 //typedef Sacado::ELRCacheFad::DFad<double> FadType;
131 //typedef Sacado::Fad::SFad<double,8> FadType;
132 typedef Sacado::CacheFad::SFad<double,8> FadType;
133 //typedef Sacado::ELRCacheFad::SFad<double,8> FadType;
134 //typedef Sacado::Fad::SLFad<double,8> FadType;
135 //typedef Sacado::CacheFad::SLFad<double,8> FadType;
136 //typedef Sacado::ELRCacheFad::SLFad<double,8> FadType;
137 
138 //#define DUMP_DATA
139 
140 // Functions to evaluate exact solution and derivatives
141 double evalu(double & x, double & y, double & z);
142 int evalGradu(double & x, double & y, double & z, double & gradu1, double & gradu2, double & gradu3);
143 double evalDivGradu(double & x, double & y, double & z);
144 
145 int main(int argc, char *argv[]) {
146 
147  // Check number of arguments
148  if (argc < 4) {
149  std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n";
150  std::cout <<"Usage:\n\n";
151  std::cout <<" ./Intrepid_example_Drivers_Example_03AD.exe NX NY NZ verbose\n\n";
152  std::cout <<" where \n";
153  std::cout <<" int NX - num intervals in x direction (assumed box domain, 0,1) \n";
154  std::cout <<" int NY - num intervals in y direction (assumed box domain, 0,1) \n";
155  std::cout <<" int NZ - num intervals in z direction (assumed box domain, 0,1) \n";
156  std::cout <<" verbose (optional) - any character, indicates verbose output \n\n";
157  exit(1);
158  }
159 
160  // This little trick lets us print to std::cout only if
161  // a (dummy) command-line argument is provided.
162  int iprint = argc - 1;
163  Teuchos::RCP<std::ostream> outStream;
164  Teuchos::oblackholestream bhs; // outputs nothing
165  if (iprint > 3)
166  outStream = Teuchos::rcp(&std::cout, false);
167  else
168  outStream = Teuchos::rcp(&bhs, false);
169 
170  // Save the format state of the original std::cout.
171  Teuchos::oblackholestream oldFormatState;
172  oldFormatState.copyfmt(std::cout);
173 
174  *outStream \
175  << "===============================================================================\n" \
176  << "| |\n" \
177  << "| Example: Generate Stiffness Matrix and Right Hand Side Vector for |\n" \
178  << "| Poisson Equation on Hexahedral Mesh |\n" \
179  << "| |\n" \
180  << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
181  << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
182  << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
183  << "| |\n" \
184  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
185  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
186  << "| |\n" \
187  << "===============================================================================\n";
188 
189 
190  // ************************************ GET INPUTS **************************************
191 
192  int NX = atoi(argv[1]); // num intervals in x direction (assumed box domain, 0,1)
193  int NY = atoi(argv[2]); // num intervals in y direction (assumed box domain, 0,1)
194  int NZ = atoi(argv[3]); // num intervals in z direction (assumed box domain, 0,1)
195 
196  // *********************************** CELL TOPOLOGY **********************************
197 
198  // Get cell topology for base hexahedron
199  typedef shards::CellTopology CellTopology;
200  CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
201 
202  // Get dimensions
203  int numNodesPerElem = hex_8.getNodeCount();
204  int spaceDim = hex_8.getDimension();
205 
206  // *********************************** GENERATE MESH ************************************
207 
208  *outStream << "Generating mesh ... \n\n";
209 
210  *outStream << " NX" << " NY" << " NZ\n";
211  *outStream << std::setw(5) << NX <<
212  std::setw(5) << NY <<
213  std::setw(5) << NZ << "\n\n";
214 
215  // Print mesh information
216  int numElems = NX*NY*NZ;
217  int numNodes = (NX+1)*(NY+1)*(NZ+1);
218  *outStream << " Number of Elements: " << numElems << " \n";
219  *outStream << " Number of Nodes: " << numNodes << " \n\n";
220 
221  // Cube
222  double leftX = 0.0, rightX = 1.0;
223  double leftY = 0.0, rightY = 1.0;
224  double leftZ = 0.0, rightZ = 1.0;
225 
226  // Mesh spacing
227  double hx = (rightX-leftX)/((double)NX);
228  double hy = (rightY-leftY)/((double)NY);
229  double hz = (rightZ-leftZ)/((double)NZ);
230 
231  // Get nodal coordinates
232  FieldContainer<double> nodeCoord(numNodes, spaceDim);
233  FieldContainer<int> nodeOnBoundary(numNodes);
234  int inode = 0;
235  for (int k=0; k<NZ+1; k++) {
236  for (int j=0; j<NY+1; j++) {
237  for (int i=0; i<NX+1; i++) {
238  nodeCoord(inode,0) = leftX + (double)i*hx;
239  nodeCoord(inode,1) = leftY + (double)j*hy;
240  nodeCoord(inode,2) = leftZ + (double)k*hz;
241  if (k==0 || j==0 || i==0 || k==NZ || j==NY || i==NX){
242  nodeOnBoundary(inode)=1;
243  }
244  else {
245  nodeOnBoundary(inode)=0;
246  }
247  inode++;
248  }
249  }
250  }
251 
252 #ifdef DUMP_DATA
253  // Print nodal coords
254  ofstream fcoordout("coords.dat");
255  for (int i=0; i<numNodes; i++) {
256  fcoordout << nodeCoord(i,0) <<" ";
257  fcoordout << nodeCoord(i,1) <<" ";
258  fcoordout << nodeCoord(i,2) <<"\n";
259  }
260  fcoordout.close();
261 #endif
262 
263 
264  // Element to Node map
265  FieldContainer<int> elemToNode(numElems, numNodesPerElem);
266  int ielem = 0;
267  for (int k=0; k<NZ; k++) {
268  for (int j=0; j<NY; j++) {
269  for (int i=0; i<NX; i++) {
270  elemToNode(ielem,0) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i;
271  elemToNode(ielem,1) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i + 1;
272  elemToNode(ielem,2) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i + 1;
273  elemToNode(ielem,3) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i;
274  elemToNode(ielem,4) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i;
275  elemToNode(ielem,5) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i + 1;
276  elemToNode(ielem,6) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i + 1;
277  elemToNode(ielem,7) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i;
278  ielem++;
279  }
280  }
281  }
282 #ifdef DUMP_DATA
283  // Output connectivity
284  ofstream fe2nout("elem2node.dat");
285  for (int k=0; k<NZ; k++) {
286  for (int j=0; j<NY; j++) {
287  for (int i=0; i<NX; i++) {
288  int ielem = i + j * NX + k * NX * NY;
289  for (int m=0; m<numNodesPerElem; m++){
290  fe2nout << elemToNode(ielem,m) <<" ";
291  }
292  fe2nout <<"\n";
293  }
294  }
295  }
296  fe2nout.close();
297 #endif
298 
299 
300  // ************************************ CUBATURE **************************************
301 
302  *outStream << "Getting cubature ... \n\n";
303 
304  // Get numerical integration points and weights
305  DefaultCubatureFactory<double> cubFactory;
306  int cubDegree = 2;
307  Teuchos::RCP<Cubature<double> > hexCub = cubFactory.create(hex_8, cubDegree);
308 
309  int cubDim = hexCub->getDimension();
310  int numCubPoints = hexCub->getNumPoints();
311 
312  FieldContainer<double> cubPoints(numCubPoints, cubDim);
313  FieldContainer<double> cubWeights(numCubPoints);
314 
315  hexCub->getCubature(cubPoints, cubWeights);
316 
317 
318  // ************************************** BASIS ***************************************
319 
320  *outStream << "Getting basis ... \n\n";
321 
322  // Define basis
324  int numFieldsG = hexHGradBasis.getCardinality();
325  FieldContainer<double> hexGVals(numFieldsG, numCubPoints);
326  FieldContainer<double> hexGrads(numFieldsG, numCubPoints, spaceDim);
327 
328  // Evaluate basis values and gradients at cubature points
329  hexHGradBasis.getValues(hexGVals, cubPoints, OPERATOR_VALUE);
330  hexHGradBasis.getValues(hexGrads, cubPoints, OPERATOR_GRAD);
331 
332 
333  // ******** FEM ASSEMBLY *************
334 
335  *outStream << "Building stiffness matrix and right hand side ... \n\n";
336 
337  // Settings and data structures for mass and stiffness matrices
339  typedef FunctionSpaceTools fst;
340  int numCells = BATCH_SIZE;
341  int numBatches = numElems/numCells;
342 
343  // Container for nodes
344  FieldContainer<double> hexNodes(numCells, numNodesPerElem, spaceDim);
345  // Containers for Jacobian
346  FieldContainer<double> hexJacobian(numCells, numCubPoints, spaceDim, spaceDim);
347  FieldContainer<double> hexJacobInv(numCells, numCubPoints, spaceDim, spaceDim);
348  FieldContainer<double> hexJacobDet(numCells, numCubPoints);
349  // Containers for element HGRAD stiffness matrix
350  FieldContainer<double> localStiffMatrix(numCells, numFieldsG, numFieldsG);
351  FieldContainer<double> weightedMeasure(numCells, numCubPoints);
352  FieldContainer<double> hexGradsTransformed(numCells, numFieldsG, numCubPoints, spaceDim);
353  FieldContainer<double> hexGradsTransformedWeighted(numCells, numFieldsG, numCubPoints, spaceDim);
354  // Containers for right hand side vectors
355  FieldContainer<double> rhsData(numCells, numCubPoints);
356  FieldContainer<double> localRHS(numCells, numFieldsG);
357  FieldContainer<double> hexGValsTransformed(numCells, numFieldsG, numCubPoints);
358  FieldContainer<double> hexGValsTransformedWeighted(numCells, numFieldsG, numCubPoints);
359  // Container for cubature points in physical space
360  FieldContainer<double> physCubPoints(numCells, numCubPoints, cubDim);
361 
362  // Global arrays in Epetra format
363  Epetra_SerialComm Comm;
364  Epetra_Map globalMapG(numNodes, 0, Comm);
365  Epetra_FECrsMatrix StiffMatrix(Copy, globalMapG, 64);
366  Epetra_FEVector rhs(globalMapG);
367  Epetra_FEVector rhsViaAD(globalMapG);
368 
369  // Additional arrays used in AD-based assembly.
370  FieldContainer<FadType> cellResidualAD(numCells, numFieldsG);
371  FieldContainer<FadType> FEFunc(numCells, numCubPoints, spaceDim);
372  FieldContainer<FadType> x_fad(numCells, numFieldsG);
373  for (int ci=0; ci<numCells; ci++) {
374  for(int j=0; j<numFieldsG; j++) {
375  x_fad(ci,j) = FadType(numFieldsG, j, 0.0);
376  }
377  }
378 
379  Teuchos::Time timer_jac_analytic("Time to compute element PDE Jacobians analytically: ");
380  Teuchos::Time timer_jac_fad ("Time to compute element PDE Jacobians using AD: ");
381  Teuchos::Time timer_jac_insert ("Time for global insert, w/o graph: ");
382  Teuchos::Time timer_jac_insert_g("Time for global insert, w/ graph: ");
383  Teuchos::Time timer_jac_ga ("Time for GlobalAssemble, w/o graph: ");
384  Teuchos::Time timer_jac_ga_g ("Time for GlobalAssemble, w/ graph: ");
385  Teuchos::Time timer_jac_fc ("Time for FillComplete, w/o graph: ");
386  Teuchos::Time timer_jac_fc_g ("Time for FillComplete, w/ graph: ");
387 
388 
389 
390 
391  // *** Analytic element loop ***
392  for (int bi=0; bi<numBatches; bi++) {
393 
394  // Physical cell coordinates
395  for (int ci=0; ci<numCells; ci++) {
396  int k = bi*numCells+ci;
397  for (int i=0; i<numNodesPerElem; i++) {
398  hexNodes(ci,i,0) = nodeCoord(elemToNode(k,i),0);
399  hexNodes(ci,i,1) = nodeCoord(elemToNode(k,i),1);
400  hexNodes(ci,i,2) = nodeCoord(elemToNode(k,i),2);
401  }
402  }
403 
404  // Compute cell Jacobians, their inverses and their determinants
405  CellTools::setJacobian(hexJacobian, cubPoints, hexNodes, hex_8);
406  CellTools::setJacobianInv(hexJacobInv, hexJacobian );
407  CellTools::setJacobianDet(hexJacobDet, hexJacobian );
408 
409  // ******************** COMPUTE ELEMENT HGrad STIFFNESS MATRICES WITHOUT AD *******************
410 
411  // transform to physical coordinates
412  fst::HGRADtransformGRAD<double>(hexGradsTransformed, hexJacobInv, hexGrads);
413 
414  // compute weighted measure
415  fst::computeCellMeasure<double>(weightedMeasure, hexJacobDet, cubWeights);
416 
417  // multiply values with weighted measure
418  fst::multiplyMeasure<double>(hexGradsTransformedWeighted,
419  weightedMeasure, hexGradsTransformed);
420 
421  // integrate to compute element stiffness matrix
422  timer_jac_analytic.start();
423  fst::integrate<double>(localStiffMatrix,
424  hexGradsTransformed, hexGradsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE);
425  timer_jac_analytic.stop();
426 
427  // assemble into global matrix
428  for (int ci=0; ci<numCells; ci++) {
429  int k = bi*numCells+ci;
430  std::vector<int> rowIndex(numFieldsG);
431  std::vector<int> colIndex(numFieldsG);
432  for (int row = 0; row < numFieldsG; row++){
433  rowIndex[row] = elemToNode(k,row);
434  }
435  for (int col = 0; col < numFieldsG; col++){
436  colIndex[col] = elemToNode(k,col);
437  }
438  // We can insert an entire matrix at a time, but we opt for rows only.
439  //timer_jac_insert.start();
440  //StiffMatrix.InsertGlobalValues(numFieldsG, &rowIndex[0], numFieldsG, &colIndex[0], &localStiffMatrix(ci,0,0));
441  //timer_jac_insert.stop();
442  for (int row = 0; row < numFieldsG; row++){
443  timer_jac_insert.start();
444  StiffMatrix.InsertGlobalValues(1, &rowIndex[row], numFieldsG, &colIndex[0], &localStiffMatrix(ci,row,0));
445  timer_jac_insert.stop();
446  }
447  }
448 
449  // ******************* COMPUTE RIGHT-HAND SIDE WITHOUT AD *******************
450 
451  // transform integration points to physical points
452  CellTools::mapToPhysicalFrame(physCubPoints, cubPoints, hexNodes, hex_8);
453 
454  // evaluate right hand side function at physical points
455  for (int ci=0; ci<numCells; ci++) {
456  for (int nPt = 0; nPt < numCubPoints; nPt++){
457  double x = physCubPoints(ci,nPt,0);
458  double y = physCubPoints(ci,nPt,1);
459  double z = physCubPoints(ci,nPt,2);
460  rhsData(ci,nPt) = evalDivGradu(x, y, z);
461  }
462  }
463 
464  // transform basis values to physical coordinates
465  fst::HGRADtransformVALUE<double>(hexGValsTransformed, hexGVals);
466 
467  // multiply values with weighted measure
468  fst::multiplyMeasure<double>(hexGValsTransformedWeighted,
469  weightedMeasure, hexGValsTransformed);
470 
471  // integrate rhs term
472  fst::integrate<double>(localRHS, rhsData, hexGValsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE);
473 
474  // assemble into global vector
475  for (int ci=0; ci<numCells; ci++) {
476  int k = bi*numCells+ci;
477  for (int row = 0; row < numFieldsG; row++){
478  int rowIndex = elemToNode(k,row);
479  double val = -localRHS(ci,row);
480  rhs.SumIntoGlobalValues(1, &rowIndex, &val);
481  }
482  }
483 
484  } // *** end analytic element loop ***
485 
486  // Assemble global objects
487  timer_jac_ga.start(); StiffMatrix.GlobalAssemble(); timer_jac_ga.stop();
488  timer_jac_fc.start(); StiffMatrix.FillComplete(); timer_jac_fc.stop();
489  rhs.GlobalAssemble();
490 
491 
492 
493 
494  // *** AD element loop ***
495 
496  Epetra_CrsGraph mgraph = StiffMatrix.Graph();
497  Epetra_FECrsMatrix StiffMatrixViaAD(Copy, mgraph);
498  //Epetra_FECrsMatrix StiffMatrixViaAD(Copy, globalMapG, numFieldsG*numFieldsG*numFieldsG);
499 
500  for (int bi=0; bi<numBatches; bi++) {
501 
502  // ******************** COMPUTE ELEMENT HGrad STIFFNESS MATRICES AND RIGHT-HAND SIDE WITH AD ********************
503 
504  // Physical cell coordinates
505  for (int ci=0; ci<numCells; ci++) {
506  int k = bi*numCells+ci;
507  for (int i=0; i<numNodesPerElem; i++) {
508  hexNodes(ci,i,0) = nodeCoord(elemToNode(k,i),0);
509  hexNodes(ci,i,1) = nodeCoord(elemToNode(k,i),1);
510  hexNodes(ci,i,2) = nodeCoord(elemToNode(k,i),2);
511  }
512  }
513 
514  // Compute cell Jacobians, their inverses and their determinants
515  CellTools::setJacobian(hexJacobian, cubPoints, hexNodes, hex_8);
516  CellTools::setJacobianInv(hexJacobInv, hexJacobian );
517  CellTools::setJacobianDet(hexJacobDet, hexJacobian );
518 
519  // transform to physical coordinates
520  fst::HGRADtransformGRAD<double>(hexGradsTransformed, hexJacobInv, hexGrads);
521 
522  // compute weighted measure
523  fst::computeCellMeasure<double>(weightedMeasure, hexJacobDet, cubWeights);
524 
525  // multiply values with weighted measure
526  fst::multiplyMeasure<double>(hexGradsTransformedWeighted, weightedMeasure, hexGradsTransformed);
527 
528  // transform integration points to physical points
529  CellTools::mapToPhysicalFrame(physCubPoints, cubPoints, hexNodes, hex_8);
530 
531  // evaluate right hand side function at physical points
532  for (int ci=0; ci<numCells; ci++) {
533  for (int nPt = 0; nPt < numCubPoints; nPt++){
534  double x = physCubPoints(ci,nPt,0);
535  double y = physCubPoints(ci,nPt,1);
536  double z = physCubPoints(ci,nPt,2);
537  rhsData(ci,nPt) = evalDivGradu(x, y, z);
538  }
539  }
540 
541  // transform basis values to physical coordinates
542  fst::HGRADtransformVALUE<double>(hexGValsTransformed, hexGVals);
543 
544  // multiply values with weighted measure
545  fst::multiplyMeasure<double>(hexGValsTransformedWeighted,
546  weightedMeasure, hexGValsTransformed);
547 
548  timer_jac_fad.start();
549  // must zero out FEFunc due to the strange default sum-into behavior of evaluate function
550  FEFunc.initialize();
551 
552  // compute FEFunc, a linear combination of the gradients of the basis functions, with coefficients x_fad
553  // this will replace the gradient of the trial function in the weak form of the equation
554  fst::evaluate<FadType>(FEFunc,x_fad,hexGradsTransformed);
555 
556  // integrate to compute element residual
557  fst::integrate<FadType>(cellResidualAD, FEFunc, hexGradsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE);
558  timer_jac_fad.stop();
559  fst::integrate<FadType>(cellResidualAD, rhsData, hexGValsTransformedWeighted, INTREPID_INTEGRATE_COMP_ENGINE, true);
560 
561  // assemble into global matrix
562  for (int ci=0; ci<numCells; ci++) {
563  int k = bi*numCells+ci;
564  std::vector<int> rowIndex(numFieldsG);
565  std::vector<int> colIndex(numFieldsG);
566  for (int row = 0; row < numFieldsG; row++){
567  rowIndex[row] = elemToNode(k,row);
568  }
569  for (int col = 0; col < numFieldsG; col++){
570  colIndex[col] = elemToNode(k,col);
571  }
572  for (int row = 0; row < numFieldsG; row++){
573  timer_jac_insert_g.start();
574  StiffMatrixViaAD.SumIntoGlobalValues(1, &rowIndex[row], numFieldsG, &colIndex[0], cellResidualAD(ci,row).dx());
575  timer_jac_insert_g.stop();
576  }
577  }
578 
579  // assemble into global vector
580  for (int ci=0; ci<numCells; ci++) {
581  int k = bi*numCells+ci;
582  for (int row = 0; row < numFieldsG; row++){
583  int rowIndex = elemToNode(k,row);
584  double val = -cellResidualAD(ci,row).val();
585  rhsViaAD.SumIntoGlobalValues(1, &rowIndex, &val);
586  }
587  }
588 
589  } // *** end AD element loop ***
590 
591  // Assemble global objects
592  timer_jac_ga_g.start(); StiffMatrixViaAD.GlobalAssemble(); timer_jac_ga_g.stop();
593  timer_jac_fc_g.start(); StiffMatrixViaAD.FillComplete(); timer_jac_fc_g.stop();
594  rhsViaAD.GlobalAssemble();
595 
596 
597 
598  /****** Output *******/
599 
600 #ifdef DUMP_DATA
601  // Dump matrices to disk
602  EpetraExt::RowMatrixToMatlabFile("stiff_matrix.dat",StiffMatrix);
603  EpetraExt::RowMatrixToMatlabFile("stiff_matrixAD.dat",StiffMatrixViaAD);
604  EpetraExt::MultiVectorToMatrixMarketFile("rhs_vector.dat",rhs,0,0,false);
605  EpetraExt::MultiVectorToMatrixMarketFile("rhs_vectorAD.dat",rhsViaAD,0,0,false);
606 #endif
607 
608  // take the infinity norm of the difference between StiffMatrix and StiffMatrixViaAD to see that
609  // the two matrices are the same
610  EpetraExt::MatrixMatrix::Add(StiffMatrix, false, 1.0, StiffMatrixViaAD, -1.0);
611  double normMat = StiffMatrixViaAD.NormInf();
612  *outStream << "Infinity norm of difference between stiffness matrices = " << normMat << "\n";
613 
614  // take the infinity norm of the difference between rhs and rhsViaAD to see that
615  // the two vectors are the same
616  double normVec;
617  rhsViaAD.Update(-1.0, rhs, 1.0);
618  rhsViaAD.NormInf(&normVec);
619  *outStream << "Infinity norm of difference between right-hand side vectors = " << normVec << "\n";
620 
621  *outStream << "\n\nNumber of global nonzeros: " << StiffMatrix.NumGlobalNonzeros() << "\n\n";
622 
623  *outStream << timer_jac_analytic.name() << " " << timer_jac_analytic.totalElapsedTime() << " sec\n";
624  *outStream << timer_jac_fad.name() << " " << timer_jac_fad.totalElapsedTime() << " sec\n\n";
625  *outStream << timer_jac_insert.name() << " " << timer_jac_insert.totalElapsedTime() << " sec\n";
626  *outStream << timer_jac_insert_g.name() << " " << timer_jac_insert_g.totalElapsedTime() << " sec\n\n";
627  *outStream << timer_jac_ga.name() << " " << timer_jac_ga.totalElapsedTime() << " sec\n";
628  *outStream << timer_jac_ga_g.name() << " " << timer_jac_ga_g.totalElapsedTime() << " sec\n\n";
629  *outStream << timer_jac_fc.name() << " " << timer_jac_fc.totalElapsedTime() << " sec\n";
630  *outStream << timer_jac_fc_g.name() << " " << timer_jac_fc_g.totalElapsedTime() << " sec\n\n";
631 
632  // Adjust stiffness matrix and rhs based on boundary conditions
633  /* skip this part ...
634  for (int row = 0; row<numNodes; row++){
635  if (nodeOnBoundary(row)) {
636  int rowindex = row;
637  for (int col=0; col<numNodes; col++){
638  double val = 0.0;
639  int colindex = col;
640  StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &colindex, &val);
641  }
642  double val = 1.0;
643  StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &rowindex, &val);
644  val = 0.0;
645  rhs.ReplaceGlobalValues(1, &rowindex, &val);
646  }
647  }
648  */
649 
650  if ((normMat < 1.0e4*INTREPID_TOL) && (normVec < 1.0e4*INTREPID_TOL)) {
651  std::cout << "End Result: TEST PASSED\n";
652  }
653  else {
654  std::cout << "End Result: TEST FAILED\n";
655  }
656 
657  // reset format state of std::cout
658  std::cout.copyfmt(oldFormatState);
659 
660  return 0;
661 }
662 
663 
664 // Calculates Laplacian of exact solution u
665  double evalDivGradu(double & x, double & y, double & z)
666  {
667  /*
668  // function 1
669  double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z);
670  */
671 
672  // function 2
673  double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
674  + 2.0*M_PI*cos(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
675  + 2.0*M_PI*cos(M_PI*y)*sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z)
676  + 2.0*M_PI*cos(M_PI*z)*sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z)
677  + 3.0*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
678 
679  return divGradu;
680  }
virtual int getCardinality() const
Returns cardinality of the basis.
Header file for the Intrepid::CellTools class.
Header file for utility class to provide multidimensional containers.
Header file for utility class to provide array tools, such as tensor contractions, etc.
Defines expert-level interfaces for the evaluation of functions and operators in physical space (supp...
Intrepid utilities.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Header file for the Intrepid::FunctionSpaceTools class.
Header file for classes providing basic linear algebra functionality in 1D, 2D and 3D...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell...
A factory class that generates specific instances of cubatures.
Teuchos::RCP< Cubature< Scalar, ArrayPoint, ArrayWeight > > create(const shards::CellTopology &cellTopology, const std::vector< int > &degree)
Factory method.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Hexahedron cell.
A stateless class for operations on cell data. Provides methods for: