85 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp" 
   91 #include "Epetra_Time.h" 
   92 #include "Epetra_Map.h" 
   93 #include "Epetra_FECrsMatrix.h" 
   94 #include "Epetra_FEVector.h" 
   95 #include "Epetra_SerialComm.h" 
   98 #include "Teuchos_oblackholestream.hpp" 
   99 #include "Teuchos_RCP.hpp" 
  100 #include "Teuchos_BLAS.hpp" 
  103 #include "Shards_CellTopology.hpp" 
  106 #include "EpetraExt_RowMatrixOut.h" 
  107 #include "EpetraExt_MultiVectorOut.h" 
  110 using namespace Intrepid;
 
  113 double evalu(
double & x, 
double & y, 
double & z);
 
  114 int evalGradu(
double & x, 
double & y, 
double & z, 
double & gradu1, 
double & gradu2, 
double & gradu3);
 
  115 double evalDivGradu(
double & x, 
double & y, 
double & z);
 
  117 int main(
int argc, 
char *argv[]) {
 
  121       std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
 
  122       std::cout <<
"Usage:\n\n";
 
  123       std::cout <<
"  ./Intrepid_example_Drivers_Example_03.exe NX NY NZ verbose\n\n";
 
  124       std::cout <<
" where \n";
 
  125       std::cout <<
"   int NX              - num intervals in x direction (assumed box domain, 0,1) \n";
 
  126       std::cout <<
"   int NY              - num intervals in y direction (assumed box domain, 0,1) \n";
 
  127       std::cout <<
"   int NZ              - num intervals in z direction (assumed box domain, 0,1) \n";
 
  128       std::cout <<
"   verbose (optional)  - any character, indicates verbose output \n\n";
 
  134   int iprint     = argc - 1;
 
  135   Teuchos::RCP<std::ostream> outStream;
 
  136   Teuchos::oblackholestream bhs; 
 
  138     outStream = Teuchos::rcp(&std::cout, 
false);
 
  140     outStream = Teuchos::rcp(&bhs, 
false);
 
  143   Teuchos::oblackholestream oldFormatState;
 
  144   oldFormatState.copyfmt(std::cout);
 
  147     << 
"===============================================================================\n" \
 
  149     << 
"|  Example: Generate Stiffness Matrix and Right Hand Side Vector for          |\n" \
 
  150     << 
"|                   Poisson Equation on Hexahedral Mesh                       |\n" \
 
  152     << 
"|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
 
  153     << 
"|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
 
  154     << 
"|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
 
  156     << 
"|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
 
  157     << 
"|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
 
  159     << 
"===============================================================================\n";
 
  164     int NX            = atoi(argv[1]);  
 
  165     int NY            = atoi(argv[2]);  
 
  166     int NZ            = atoi(argv[3]);  
 
  171     typedef shards::CellTopology    CellTopology;
 
  172     CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
 
  175     int numNodesPerElem = hex_8.getNodeCount();
 
  176     int spaceDim = hex_8.getDimension();
 
  180     *outStream << 
"Generating mesh ... \n\n";
 
  182     *outStream << 
"   NX" << 
"   NY" << 
"   NZ\n";
 
  183     *outStream << std::setw(5) << NX <<
 
  184                  std::setw(5) << NY <<
 
  185                  std::setw(5) << NZ << 
"\n\n";
 
  188     int numElems = NX*NY*NZ;
 
  189     int numNodes = (NX+1)*(NY+1)*(NZ+1);
 
  190     *outStream << 
" Number of Elements: " << numElems << 
" \n";
 
  191     *outStream << 
"    Number of Nodes: " << numNodes << 
" \n\n";
 
  194     double leftX = 0.0, rightX = 1.0;
 
  195     double leftY = 0.0, rightY = 1.0;
 
  196     double leftZ = 0.0, rightZ = 1.0;
 
  199     double hx = (rightX-leftX)/((
double)NX);
 
  200     double hy = (rightY-leftY)/((
double)NY);
 
  201     double hz = (rightZ-leftZ)/((
double)NZ);
 
  207     for (
int k=0; k<NZ+1; k++) {
 
  208       for (
int j=0; j<NY+1; j++) {
 
  209         for (
int i=0; i<NX+1; i++) {
 
  210           nodeCoord(inode,0) = leftX + (double)i*hx;
 
  211           nodeCoord(inode,1) = leftY + (double)j*hy;
 
  212           nodeCoord(inode,2) = leftZ + (double)k*hz;
 
  213           if (k==0 || j==0 || i==0 || k==NZ || j==NY || i==NX){
 
  214              nodeOnBoundary(inode)=1;
 
  217              nodeOnBoundary(inode)=0;
 
  226     ofstream fcoordout(
"coords.dat");
 
  227     for (
int i=0; i<numNodes; i++) {
 
  228        fcoordout << nodeCoord(i,0) <<
" ";
 
  229        fcoordout << nodeCoord(i,1) <<
" ";
 
  230        fcoordout << nodeCoord(i,2) <<
"\n";
 
  239     for (
int k=0; k<NZ; k++) {
 
  240       for (
int j=0; j<NY; j++) {
 
  241         for (
int i=0; i<NX; i++) {
 
  242           elemToNode(ielem,0) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i;
 
  243           elemToNode(ielem,1) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i + 1;
 
  244           elemToNode(ielem,2) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i + 1;
 
  245           elemToNode(ielem,3) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i;
 
  246           elemToNode(ielem,4) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i;
 
  247           elemToNode(ielem,5) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i + 1;
 
  248           elemToNode(ielem,6) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i + 1;
 
  249           elemToNode(ielem,7) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i;
 
  256     ofstream fe2nout(
"elem2node.dat");
 
  257     for (
int k=0; k<NZ; k++) {
 
  258       for (
int j=0; j<NY; j++) {
 
  259         for (
int i=0; i<NX; i++) {
 
  260           ielem = i + j * NX + k * NX * NY;
 
  261           for (
int m=0; m<numNodesPerElem; m++){
 
  262               fe2nout << elemToNode(ielem,m) <<
"  ";
 
  274     *outStream << 
"Getting cubature ... \n\n";
 
  279     Teuchos::RCP<Cubature<double> > hexCub = cubFactory.
create(hex_8, cubDegree); 
 
  281     int cubDim       = hexCub->getDimension();
 
  282     int numCubPoints = hexCub->getNumPoints();
 
  287     hexCub->getCubature(cubPoints, cubWeights);
 
  292      *outStream << 
"Getting basis ... \n\n";
 
  301      hexHGradBasis.
getValues(hexGVals, cubPoints, OPERATOR_VALUE);
 
  302      hexHGradBasis.
getValues(hexGrads, cubPoints, OPERATOR_GRAD);
 
  307     *outStream << 
"Building stiffness matrix and right hand side ... \n\n";
 
  334     Epetra_SerialComm Comm;
 
  335     Epetra_Map globalMapG(numNodes, 0, Comm);
 
  336     Epetra_FECrsMatrix StiffMatrix(Copy, globalMapG, numFieldsG);
 
  337     Epetra_FEVector rhs(globalMapG);
 
  340     for (
int k=0; k<numElems; k++) {
 
  343       for (
int i=0; i<numNodesPerElem; i++) {
 
  344          hexNodes(0,i,0) = nodeCoord(elemToNode(k,i),0);
 
  345          hexNodes(0,i,1) = nodeCoord(elemToNode(k,i),1);
 
  346          hexNodes(0,i,2) = nodeCoord(elemToNode(k,i),2);
 
  350        CellTools::setJacobian(hexJacobian, cubPoints, hexNodes, hex_8);
 
  351        CellTools::setJacobianInv(hexJacobInv, hexJacobian );
 
  352        CellTools::setJacobianDet(hexJacobDet, hexJacobian );
 
  357       fst::HGRADtransformGRAD<double>(hexGradsTransformed, hexJacobInv, hexGrads);
 
  360       fst::computeCellMeasure<double>(weightedMeasure, hexJacobDet, cubWeights);
 
  363       fst::multiplyMeasure<double>(hexGradsTransformedWeighted,
 
  364                                    weightedMeasure, hexGradsTransformed);
 
  367       fst::integrate<double>(localStiffMatrix,
 
  368                              hexGradsTransformed, hexGradsTransformedWeighted, COMP_BLAS);
 
  371       for (
int row = 0; row < numFieldsG; row++){
 
  372         for (
int col = 0; col < numFieldsG; col++){
 
  373             int rowIndex = elemToNode(k,row);
 
  374             int colIndex = elemToNode(k,col);
 
  375             double val = localStiffMatrix(0,row,col);
 
  376             StiffMatrix.InsertGlobalValues(1, &rowIndex, 1, &colIndex, &val);
 
  383        CellTools::mapToPhysicalFrame(physCubPoints, cubPoints, hexNodes, hex_8);
 
  386        for (
int nPt = 0; nPt < numCubPoints; nPt++){
 
  388           double x = physCubPoints(0,nPt,0);
 
  389           double y = physCubPoints(0,nPt,1);
 
  390           double z = physCubPoints(0,nPt,2);
 
  392           rhsData(0,nPt) = evalDivGradu(x, y, z);
 
  396       fst::HGRADtransformVALUE<double>(hexGValsTransformed, hexGVals);
 
  399       fst::multiplyMeasure<double>(hexGValsTransformedWeighted,
 
  400                                    weightedMeasure, hexGValsTransformed);
 
  403       fst::integrate<double>(localRHS, rhsData, hexGValsTransformedWeighted, 
 
  407      for (
int row = 0; row < numFieldsG; row++){
 
  408            int rowIndex = elemToNode(k,row);
 
  409            double val = -localRHS(0,row);
 
  410            rhs.SumIntoGlobalValues(1, &rowIndex, &val);
 
  417    StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
 
  418    rhs.GlobalAssemble();
 
  422    for (
int row = 0; row<numNodes; row++){
 
  423        if (nodeOnBoundary(row)) {
 
  425           for (
int col=0; col<numNodes; col++){
 
  428               StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &colindex, &val);
 
  431           StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &rowindex, &val);
 
  433           rhs.ReplaceGlobalValues(1, &rowindex, &val);
 
  439      EpetraExt::RowMatrixToMatlabFile(
"stiff_matrix.dat",StiffMatrix);
 
  440      EpetraExt::MultiVectorToMatrixMarketFile(
"rhs_vector.dat",rhs,0,0,
false);
 
  443    std::cout << 
"End Result: TEST PASSED\n";
 
  446    std::cout.copyfmt(oldFormatState);
 
  453  double evalu(
double & x, 
double & y, 
double & z)
 
  461    double exactu = sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
 
  467  int evalGradu(
double & x, 
double & y, 
double & z, 
double & gradu1, 
double & gradu2, 
double & gradu3)
 
  477        gradu1 = (M_PI*cos(M_PI*x)+sin(M_PI*x))
 
  478                   *sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
 
  479        gradu2 = (M_PI*cos(M_PI*y)+sin(M_PI*y))
 
  480                   *sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z);
 
  481        gradu3 = (M_PI*cos(M_PI*z)+sin(M_PI*z))
 
  482                   *sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z);
 
  488  double evalDivGradu(
double & x, 
double & y, 
double & z)
 
  496    double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
 
  497                     + 2.0*M_PI*cos(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
 
  498                     + 2.0*M_PI*cos(M_PI*y)*sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z)
 
  499                     + 2.0*M_PI*cos(M_PI*z)*sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z)
 
  500                     + 3.0*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
 
virtual int getCardinality() const 
Returns cardinality of the basis. 
Header file for utility class to provide multidimensional containers. 
Header file for the abstract base class Intrepid::DefaultCubatureFactory. 
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 > °ree)
Factory method. 
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const 
Evaluation of a FEM basis on a reference Hexahedron cell.