98 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp" 
  104 #include "Epetra_Time.h" 
  105 #include "Epetra_Map.h" 
  106 #include "Epetra_SerialComm.h" 
  107 #include "Epetra_FECrsMatrix.h" 
  108 #include "Epetra_FEVector.h" 
  109 #include "Epetra_Vector.h" 
  112 #include "Teuchos_oblackholestream.hpp" 
  113 #include "Teuchos_RCP.hpp" 
  114 #include "Teuchos_BLAS.hpp" 
  117 #include "Shards_CellTopology.hpp" 
  120 #include "EpetraExt_RowMatrixOut.h" 
  121 #include "EpetraExt_MultiVectorOut.h" 
  124 using namespace Intrepid;
 
  127 int evalu(
double & uExact0, 
double & uExact1, 
double & uExact2, 
double & x, 
double & y, 
double & z);
 
  128 double evalDivu(
double & x, 
double & y, 
double & z);
 
  129 int evalCurlu(
double & curlu0, 
double & curlu1, 
double & curlu2, 
double & x, 
double & y, 
double & z);
 
  130 int evalGradDivu(
double & gradDivu0, 
double & gradDivu1, 
double & gradDivu2, 
double & x, 
double & y, 
double & z);
 
  132 int main(
int argc, 
char *argv[]) {
 
  136       std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
 
  137       std::cout <<
"Usage:\n\n";
 
  138       std::cout <<
"  ./Intrepid_example_Drivers_Example_01.exe NX NY NZ randomMesh mu1 mu2 mu1LX mu1RX mu1LY mu1RY mu1LZ mu1RZ verbose\n\n";
 
  139       std::cout <<
" where \n";
 
  140       std::cout <<
"   int NX              - num intervals in x direction (assumed box domain, -1,1) \n";
 
  141       std::cout <<
"   int NY              - num intervals in y direction (assumed box domain, -1,1) \n";
 
  142       std::cout <<
"   int NZ              - num intervals in z direction (assumed box domain, -1,1) \n";
 
  143       std::cout <<
"   int randomMesh      - 1 if mesh randomizer is to be used 0 if not \n";
 
  144       std::cout <<
"   double mu1          - material property value for region 1 \n";
 
  145       std::cout <<
"   double mu2          - material property value for region 2 \n";
 
  146       std::cout <<
"   double mu1LX        - left X boundary for region 1 \n";
 
  147       std::cout <<
"   double mu1RX        - right X boundary for region 1 \n";
 
  148       std::cout <<
"   double mu1LY        - left Y boundary for region 1 \n";
 
  149       std::cout <<
"   double mu1RY        - right Y boundary for region 1 \n";
 
  150       std::cout <<
"   double mu1LZ        - bottom Z boundary for region 1 \n";
 
  151       std::cout <<
"   double mu1RZ        - top Z boundary for region 1 \n";
 
  152       std::cout <<
"   verbose (optional)  - any character, indicates verbose output \n\n";
 
  158   int iprint     = argc - 1;
 
  159   Teuchos::RCP<std::ostream> outStream;
 
  160   Teuchos::oblackholestream bhs; 
 
  162     outStream = Teuchos::rcp(&std::cout, 
false);
 
  164     outStream = Teuchos::rcp(&bhs, 
false);
 
  167   Teuchos::oblackholestream oldFormatState;
 
  168   oldFormatState.copyfmt(std::cout);
 
  171     << 
"===============================================================================\n" \
 
  173     << 
"|  Example: Generate Mass and Stiffness Matrices and Right-Hand Side Vector   |\n" 
  174     << 
"|    for Div-Curl System on Hexahedral Mesh with Curl-Conforming Elements     |\n" \
 
  176     << 
"|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
 
  177     << 
"|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
 
  178     << 
"|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
 
  180     << 
"|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
 
  181     << 
"|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
 
  183     << 
"===============================================================================\n";
 
  192     int NX            = atoi(argv[1]);  
 
  193     int NY            = atoi(argv[2]);  
 
  194     int NZ            = atoi(argv[3]);  
 
  195     int randomMesh    = atoi(argv[4]);  
 
  196     double mu1        = atof(argv[5]);  
 
  197     double mu2        = atof(argv[6]);  
 
  198     double mu1LeftX   = atof(argv[7]);  
 
  199     double mu1RightX  = atof(argv[8]);  
 
  200     double mu1LeftY   = atof(argv[9]);  
 
  201     double mu1RightY  = atof(argv[10]); 
 
  202     double mu1LeftZ   = atof(argv[11]); 
 
  203     double mu1RightZ  = atof(argv[12]); 
 
  208     typedef shards::CellTopology    CellTopology;
 
  209     CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
 
  212     int numNodesPerElem = hex_8.getNodeCount();
 
  213     int numEdgesPerElem = hex_8.getEdgeCount();
 
  214     int numFacesPerElem = hex_8.getSideCount();
 
  215     int numNodesPerEdge = 2;
 
  216     int numNodesPerFace = 4;
 
  217     int numEdgesPerFace = 4;
 
  218     int spaceDim = hex_8.getDimension();
 
  222     for (
int i=0; i<numEdgesPerElem; i++){
 
  223         refEdgeToNode(i,0)=hex_8.getNodeMap(1, i, 0);
 
  224         refEdgeToNode(i,1)=hex_8.getNodeMap(1, i, 1);
 
  229     *outStream << 
"Generating mesh ... \n\n";
 
  231     *outStream << 
"    NX" << 
"   NY" << 
"   NZ\n";
 
  232     *outStream << std::setw(5) << NX <<
 
  233                  std::setw(5) << NY <<
 
  234                  std::setw(5) << NZ << 
"\n\n";
 
  237     int numElems = NX*NY*NZ;
 
  238     int numNodes = (NX+1)*(NY+1)*(NZ+1);
 
  239     int numEdges = (NX)*(NY + 1)*(NZ + 1) + (NX + 1)*(NY)*(NZ + 1) + (NX + 1)*(NY + 1)*(NZ);
 
  240     int numFaces = (NX)*(NY)*(NZ + 1) + (NX)*(NY + 1)*(NZ) + (NX + 1)*(NY)*(NZ);
 
  241     *outStream << 
" Number of Elements: " << numElems << 
" \n";
 
  242     *outStream << 
"    Number of Nodes: " << numNodes << 
" \n";
 
  243     *outStream << 
"    Number of Edges: " << numEdges << 
" \n";
 
  244     *outStream << 
"    Number of Faces: " << numFaces << 
" \n\n";
 
  247     double leftX = -1.0, rightX = 1.0;
 
  248     double leftY = -1.0, rightY = 1.0;
 
  249     double leftZ = -1.0, rightZ = 1.0;
 
  252     double hx = (rightX-leftX)/((
double)NX);
 
  253     double hy = (rightY-leftY)/((
double)NY);
 
  254     double hz = (rightZ-leftZ)/((
double)NZ);
 
  260     for (
int k=0; k<NZ+1; k++) {
 
  261       for (
int j=0; j<NY+1; j++) {
 
  262         for (
int i=0; i<NX+1; i++) {
 
  263           nodeCoord(inode,0) = leftX + (double)i*hx;
 
  264           nodeCoord(inode,1) = leftY + (double)j*hy;
 
  265           nodeCoord(inode,2) = leftZ + (double)k*hz;
 
  266           if (k==0 || j==0 || i==0 || k==NZ || j==NY || i==NX){
 
  267              nodeOnBoundary(inode)=1; 
 
  277     for (
int k=0; k<NZ; k++) {
 
  278       for (
int j=0; j<NY; j++) {
 
  279         for (
int i=0; i<NX; i++) {
 
  280           elemToNode(ielem,0) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i;
 
  281           elemToNode(ielem,1) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i + 1;
 
  282           elemToNode(ielem,2) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i + 1;
 
  283           elemToNode(ielem,3) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i;
 
  284           elemToNode(ielem,4) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i;
 
  285           elemToNode(ielem,5) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i + 1;
 
  286           elemToNode(ielem,6) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i + 1;
 
  287           elemToNode(ielem,7) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i;
 
  298     for (
int k=0; k<NZ+1; k++) {
 
  299       for (
int j=0; j<NY+1; j++) {
 
  300         for (
int i=0; i<NX+1; i++) {
 
  302                edgeToNode(iedge,0) = inode;
 
  303                edgeToNode(iedge,1) = inode + 1;
 
  304                if (j < NY && k < NZ){
 
  305                   ielem=i+j*NX+k*NX*NY;
 
  306                   elemToEdge(ielem,0) = iedge;
 
  308                      elemToEdge(ielem-NX,2) = iedge; 
 
  310                      elemToEdge(ielem-NX*NY,4) = iedge; 
 
  312                      elemToEdge(ielem-NX*NY-NX,6) = iedge; 
 
  314                else if (j == NY && k == NZ){
 
  315                   ielem=i+(NY-1)*NX+(NZ-1)*NX*NY;
 
  316                   elemToEdge(ielem,6) = iedge;
 
  318                else if (k == NZ && j < NY){
 
  319                   ielem=i+j*NX+(NZ-1)*NX*NY;
 
  320                   elemToEdge(ielem,4) = iedge;
 
  322                     elemToEdge(ielem-NX,6) = iedge;
 
  324                else if (k < NZ && j == NY){
 
  325                   ielem=i+(NY-1)*NX+k*NX*NY;
 
  326                   elemToEdge(ielem,2) = iedge;
 
  328                      elemToEdge(ielem-NX*NY,6) = iedge;
 
  333                edgeToNode(iedge,0) = inode;
 
  334                edgeToNode(iedge,1) = inode + NX+1;
 
  335                if (i < NX && k < NZ){
 
  336                   ielem=i+j*NX+k*NX*NY;
 
  337                   elemToEdge(ielem,3) = iedge;
 
  339                      elemToEdge(ielem-1,1) = iedge; 
 
  341                      elemToEdge(ielem-NX*NY,7) = iedge; 
 
  343                      elemToEdge(ielem-NX*NY-1,5) = iedge; 
 
  345                else if (i == NX && k == NZ){
 
  346                   ielem=NX-1+j*NX+(NZ-1)*NX*NY;
 
  347                   elemToEdge(ielem,5) = iedge;
 
  349                else if (k == NZ && i < NX){
 
  350                   ielem=i+j*NX+(NZ-1)*NX*NY;
 
  351                   elemToEdge(ielem,7) = iedge;
 
  353                     elemToEdge(ielem-1,5) = iedge;
 
  355                else if (k < NZ && i == NX){
 
  356                   ielem=NX-1+j*NX+k*NX*NY;
 
  357                   elemToEdge(ielem,1) = iedge;
 
  359                      elemToEdge(ielem-NX*NY,5) = iedge;
 
  364                edgeToNode(iedge,0) = inode;
 
  365                edgeToNode(iedge,1) = inode + (NX+1)*(NY+1);
 
  366                if (i < NX && j < NY){
 
  367                   ielem=i+j*NX+k*NX*NY;
 
  368                   elemToEdge(ielem,8) = iedge;
 
  370                      elemToEdge(ielem-1,9) = iedge; 
 
  372                      elemToEdge(ielem-NX,11) = iedge; 
 
  374                      elemToEdge(ielem-NX-1,10) = iedge; 
 
  376                else if (i == NX && j == NY){
 
  377                   ielem=NX-1+(NY-1)*NX+k*NX*NY;
 
  378                   elemToEdge(ielem,10) = iedge;
 
  380                else if (j == NY && i < NX){
 
  381                   ielem=i+(NY-1)*NX+k*NX*NY;
 
  382                   elemToEdge(ielem,11) = iedge;
 
  384                     elemToEdge(ielem-1,10) = iedge;
 
  386                else if (j < NY && i == NX){
 
  387                   ielem=NX-1+j*NX+k*NX*NY;
 
  388                   elemToEdge(ielem,9) = iedge;
 
  390                      elemToEdge(ielem-NX,10) = iedge;
 
  401     for (
int i=0; i<numEdges; i++){
 
  402        if (nodeOnBoundary(edgeToNode(i,0)) && nodeOnBoundary(edgeToNode(i,1))){
 
  413     for (
int k=0; k<NZ+1; k++) {
 
  414       for (
int j=0; j<NY+1; j++) {
 
  415         for (
int i=0; i<NX+1; i++) {
 
  416            if (i < NX && k < NZ) {
 
  417               faceToNode(iface,0)=inode;
 
  418               faceToNode(iface,1)=inode + 1;
 
  419               faceToNode(iface,2)=inode + (NX+1)*(NY+1)+1;
 
  420               faceToNode(iface,3)=inode + (NX+1)*(NY+1);
 
  422                  ielem=i+j*NX+k*NX*NY;
 
  423                  faceToEdge(iface,0)=elemToEdge(ielem,0);
 
  424                  faceToEdge(iface,1)=elemToEdge(ielem,9);
 
  425                  faceToEdge(iface,2)=elemToEdge(ielem,4);
 
  426                  faceToEdge(iface,3)=elemToEdge(ielem,8);
 
  427                  elemToFace(ielem,0)=iface;
 
  429                     elemToFace(ielem-NX,2)=iface;
 
  433                  ielem=i+(NY-1)*NX+k*NX*NY;
 
  434                  faceToEdge(iface,0)=elemToEdge(ielem,2);
 
  435                  faceToEdge(iface,1)=elemToEdge(ielem,10);
 
  436                  faceToEdge(iface,2)=elemToEdge(ielem,6);
 
  437                  faceToEdge(iface,3)=elemToEdge(ielem,11);
 
  438                  elemToFace(ielem,2)=iface;
 
  442            if (j < NY && k < NZ) {
 
  443               faceToNode(iface,0)=inode;
 
  444               faceToNode(iface,1)=inode + NX+1;
 
  445               faceToNode(iface,2)=inode + (NX+1)*(NY+1) + NX+1;
 
  446               faceToNode(iface,3)=inode + (NX+1)*(NY+1);
 
  448                  ielem=i+j*NX+k*NX*NY;
 
  449                  faceToEdge(iface,0)=elemToEdge(ielem,3);
 
  450                  faceToEdge(iface,1)=elemToEdge(ielem,11);
 
  451                  faceToEdge(iface,2)=elemToEdge(ielem,7);
 
  452                  faceToEdge(iface,3)=elemToEdge(ielem,8);
 
  453                  elemToFace(ielem,3)=iface;
 
  455                     elemToFace(ielem-1,1)=iface;
 
  459                  ielem=NX-1+j*NX+k*NX*NY;
 
  460                  faceToEdge(iface,0)=elemToEdge(ielem,1);
 
  461                  faceToEdge(iface,1)=elemToEdge(ielem,10);
 
  462                  faceToEdge(iface,2)=elemToEdge(ielem,5);
 
  463                  faceToEdge(iface,3)=elemToEdge(ielem,9);
 
  464                  elemToFace(ielem,1)=iface;
 
  468            if (i < NX && j < NY) {
 
  469               faceToNode(iface,0)=inode;
 
  470               faceToNode(iface,1)=inode + 1;
 
  471               faceToNode(iface,2)=inode + NX+2;
 
  472               faceToNode(iface,3)=inode + NX+1;
 
  474                  ielem=i+j*NX+k*NX*NY;
 
  475                  faceToEdge(iface,0)=elemToEdge(ielem,0);
 
  476                  faceToEdge(iface,1)=elemToEdge(ielem,1);
 
  477                  faceToEdge(iface,2)=elemToEdge(ielem,2);
 
  478                  faceToEdge(iface,3)=elemToEdge(ielem,3);
 
  479                  elemToFace(ielem,4)=iface;
 
  481                     elemToFace(ielem-NX*NY,5)=iface;
 
  485                  ielem=i+j*NX+(NZ-1)*NX*NY;
 
  486                  faceToEdge(iface,0)=elemToEdge(ielem,4);
 
  487                  faceToEdge(iface,1)=elemToEdge(ielem,5);
 
  488                  faceToEdge(iface,2)=elemToEdge(ielem,6);
 
  489                  faceToEdge(iface,3)=elemToEdge(ielem,7);
 
  490                  elemToFace(ielem,5)=iface;
 
  501     for (
int i=0; i<numFaces; i++){
 
  502        if (nodeOnBoundary(faceToNode(i,0)) && nodeOnBoundary(faceToNode(i,1))
 
  503           && nodeOnBoundary(faceToNode(i,2)) && nodeOnBoundary(faceToNode(i,3))){
 
  512     ofstream fe2nout(
"elem2node.dat");
 
  513     ofstream fe2eout(
"elem2edge.dat");
 
  514     for (
int k=0; k<NZ; k++) {
 
  515       for (
int j=0; j<NY; j++) {
 
  516         for (
int i=0; i<NX; i++) {
 
  517           int ielem = i + j * NX + k * NX * NY;
 
  518           for (
int m=0; m<numNodesPerElem; m++){
 
  519               fe2nout << elemToNode(ielem,m) <<
"  ";
 
  522           for (
int l=0; l<numEdgesPerElem; l++) {
 
  523              fe2eout << elemToEdge(ielem,l) << 
"  ";
 
  533 #ifdef DUMP_DATA_EXTRA 
  534     ofstream fed2nout(
"edge2node.dat");
 
  535     for (
int i=0; i<numEdges; i++) {
 
  536        fed2nout << edgeToNode(i,0) <<
" ";
 
  537        fed2nout << edgeToNode(i,1) <<
"\n";
 
  541     ofstream fBnodeout(
"nodeOnBndy.dat");
 
  542     ofstream fBedgeout(
"edgeOnBndy.dat");
 
  543     for (
int i=0; i<numNodes; i++) {
 
  544         fBnodeout << nodeOnBoundary(i) <<
"\n";
 
  546     for (
int i=0; i<numEdges; i++) {
 
  547         fBedgeout << edgeOnBoundary(i) <<
"\n";
 
  555     for (
int k=0; k<NZ; k++) {
 
  556       for (
int j=0; j<NY; j++) {
 
  557         for (
int i=0; i<NX; i++) {
 
  558           int ielem = i + j * NX + k * NX * NY;
 
  559           double midElemX = nodeCoord(elemToNode(ielem,0),0) + hx/2.0;
 
  560           double midElemY = nodeCoord(elemToNode(ielem,0),1) + hy/2.0;
 
  561           double midElemZ = nodeCoord(elemToNode(ielem,0),2) + hz/2.0;
 
  562           if ( (midElemX > mu1LeftX) && (midElemY > mu1LeftY) && (midElemZ > mu1LeftZ) &&
 
  563                (midElemX <= mu1RightX) && (midElemY <= mu1RightY) && (midElemZ <= mu1RightZ) ){
 
  575       for (
int k=1; k<NZ; k++) {
 
  576         for (
int j=1; j<NY; j++) {
 
  577           for (
int i=1; i<NX; i++) {
 
  578             int inode = i + j * (NX + 1) + k * (NX + 1) * (NY + 1);
 
  580             double rx = 2.0 * (double)rand()/RAND_MAX - 1.0;
 
  581             double ry = 2.0 * (double)rand()/RAND_MAX - 1.0;
 
  582             double rz = 2.0 * (double)rand()/RAND_MAX - 1.0; 
 
  584             nodeCoord(inode,0) = nodeCoord(inode,0) + 0.125 * hx * rx;
 
  585             nodeCoord(inode,1) = nodeCoord(inode,1) + 0.125 * hy * ry;
 
  586             nodeCoord(inode,2) = nodeCoord(inode,2) + 0.125 * hz * rz;
 
  594     ofstream fcoordout(
"coords.dat");
 
  595     for (
int i=0; i<numNodes; i++) {
 
  596        fcoordout << nodeCoord(i,0) <<
" ";
 
  597        fcoordout << nodeCoord(i,1) <<
" ";
 
  598        fcoordout << nodeCoord(i,2) <<
"\n";
 
  607     *outStream << 
"Building incidence matrix ... \n\n";
 
  609     Epetra_SerialComm Comm;
 
  610     Epetra_Map globalMapC(numEdges, 0, Comm);
 
  611     Epetra_Map globalMapG(numNodes, 0, Comm);
 
  612     Epetra_FECrsMatrix DGrad(Copy, globalMapC, globalMapG, 2);
 
  615     vals[0]=-1.0; vals[1]=1.0;
 
  616     for (
int j=0; j<numEdges; j++){
 
  619         colNum[0] = edgeToNode(j,0);
 
  620         colNum[1] = edgeToNode(j,1);
 
  621         DGrad.InsertGlobalValues(1, &rowNum, 2, colNum, vals);
 
  628     *outStream << 
"Getting cubature ... \n\n";
 
  632     Teuchos::RCP<Cubature<double> > hexCub = cubFactory.
create(hex_8, cubDegree); 
 
  634     int cubDim       = hexCub->getDimension();
 
  635     int numCubPoints = hexCub->getNumPoints();
 
  640     hexCub->getCubature(cubPoints, cubWeights);
 
  647     CellTopology paramQuadFace(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
 
  651     Teuchos::RCP<Cubature<double> > hexFaceCubature = cubFactoryFace.
create(paramQuadFace, 3);
 
  652     int cubFaceDim    = hexFaceCubature -> getDimension();
 
  653     int numFacePoints = hexFaceCubature -> getNumPoints();
 
  660     hexFaceCubature -> getCubature(paramGaussPoints, paramGaussWeights);
 
  666     *outStream << 
"Getting basis ... \n\n";
 
  680      hexHCurlBasis.
getValues(hexCVals, cubPoints, OPERATOR_VALUE);
 
  681      hexHCurlBasis.
getValues(hexCurls, cubPoints, OPERATOR_CURL);
 
  682      hexHGradBasis.
getValues(hexGVals, cubPoints, OPERATOR_VALUE);
 
  688     *outStream << 
"Building mass and stiffness matrices ... \n\n";
 
  739     Epetra_FECrsMatrix MassG(Copy, globalMapG, numFieldsG);
 
  740     Epetra_FECrsMatrix MassC(Copy, globalMapC, numFieldsC);
 
  741     Epetra_FECrsMatrix StiffC(Copy, globalMapC, numFieldsC);
 
  742     Epetra_FEVector rhsC(globalMapC);
 
  745     ofstream fSignsout(
"edgeSigns.dat");
 
  749     for (
int k=0; k<numElems; k++) {
 
  752       for (
int i=0; i<numNodesPerElem; i++) {
 
  753          hexNodes(0,i,0) = nodeCoord(elemToNode(k,i),0);
 
  754          hexNodes(0,i,1) = nodeCoord(elemToNode(k,i),1);
 
  755          hexNodes(0,i,2) = nodeCoord(elemToNode(k,i),2);
 
  759       for (
int j=0; j<numEdgesPerElem; j++) {
 
  760           if (elemToNode(k,refEdgeToNode(j,0))==edgeToNode(elemToEdge(k,j),0) &&
 
  761               elemToNode(k,refEdgeToNode(j,1))==edgeToNode(elemToEdge(k,j),1))
 
  762               hexEdgeSigns(0,j) = 1.0;
 
  764               hexEdgeSigns(0,j) = -1.0;
 
  766           fSignsout << hexEdgeSigns(0,j) << 
"  ";
 
  774        CellTools::setJacobian(hexJacobian, cubPoints, hexNodes, hex_8);
 
  775        CellTools::setJacobianInv(hexJacobInv, hexJacobian );
 
  776        CellTools::setJacobianDet(hexJacobDet, hexJacobian );
 
  781       fst::HGRADtransformVALUE<double>(hexGValsTransformed, hexGVals);
 
  784       fst::computeCellMeasure<double>(weightedMeasure, hexJacobDet, cubWeights);
 
  787       for (
int nC = 0; nC < numCells; nC++){
 
  788         for (
int nPt = 0; nPt < numCubPoints; nPt++){
 
  789           weightedMeasureMuInv(nC,nPt) = weightedMeasure(nC,nPt) / muVal(k);
 
  794       fst::multiplyMeasure<double>(hexGValsTransformedWeighted,
 
  795                                    weightedMeasureMuInv, hexGValsTransformed);
 
  798       fst::integrate<double>(massMatrixG,
 
  799                              hexGValsTransformed, hexGValsTransformedWeighted, COMP_CPP);
 
  803       for (
int row = 0; row < numFieldsG; row++){
 
  804         for (
int col = 0; col < numFieldsG; col++){
 
  805             int rowIndex = elemToNode(k,row);
 
  806             int colIndex = elemToNode(k,col);
 
  807             double val = massMatrixG(0,row,col);
 
  808             MassG.InsertGlobalValues(1, &rowIndex, 1, &colIndex, &val);
 
  815       fst::HCURLtransformVALUE<double>(hexCValsTransformed, hexJacobInv, 
 
  819       fst::multiplyMeasure<double>(hexCValsTransformedWeighted,
 
  820                                    weightedMeasure, hexCValsTransformed);
 
  823       fst::integrate<double>(massMatrixC,
 
  824                              hexCValsTransformed, hexCValsTransformedWeighted,
 
  828       fst::applyLeftFieldSigns<double>(massMatrixC, hexEdgeSigns);
 
  829       fst::applyRightFieldSigns<double>(massMatrixC, hexEdgeSigns);
 
  834       for (
int row = 0; row < numFieldsC; row++){
 
  835         for (
int col = 0; col < numFieldsC; col++){
 
  836             int rowIndex = elemToEdge(k,row);
 
  837             int colIndex = elemToEdge(k,col);
 
  838             double val = massMatrixC(0,row,col);
 
  839             MassC.InsertGlobalValues(1, &rowIndex, 1, &colIndex, &val);
 
  846       fst::HCURLtransformCURL<double>(hexCurlsTransformed, hexJacobian, hexJacobDet, 
 
  850       for (
int nC = 0; nC < numCells; nC++){
 
  851         for (
int nPt = 0; nPt < numCubPoints; nPt++){
 
  852           weightedMeasureMu(nC,nPt) = weightedMeasure(nC,nPt) / muVal(k);
 
  857       fst::multiplyMeasure<double>(hexCurlsTransformedWeighted,
 
  858                                    weightedMeasureMu, hexCurlsTransformed);
 
  861       fst::integrate<double>(stiffMatrixC,
 
  862                              hexCurlsTransformed, hexCurlsTransformedWeighted,
 
  866      fst::applyLeftFieldSigns<double>(stiffMatrixC, hexEdgeSigns);
 
  867      fst::applyRightFieldSigns<double>(stiffMatrixC, hexEdgeSigns);
 
  871       for (
int row = 0; row < numFieldsC; row++){
 
  872         for (
int col = 0; col < numFieldsC; col++){
 
  873             int rowIndex = elemToEdge(k,row);
 
  874             int colIndex = elemToEdge(k,col);
 
  875             double val = stiffMatrixC(0,row,col);
 
  876             StiffC.InsertGlobalValues(1, &rowIndex, 1, &colIndex, &val);
 
  884        CellTools::mapToPhysicalFrame(physCubPoints, cubPoints, hexNodes, hex_8);
 
  889        for (
int nPt = 0; nPt < numCubPoints; nPt++){
 
  891           double x = physCubPoints(0,nPt,0);
 
  892           double y = physCubPoints(0,nPt,1);
 
  893           double z = physCubPoints(0,nPt,2);
 
  894           double du1, du2, du3;
 
  896           evalCurlu(du1, du2, du3, x, y, z);
 
  897           rhsDatag(0,nPt,0) = du1;
 
  898           rhsDatag(0,nPt,1) = du2;
 
  899           rhsDatag(0,nPt,2) = du3;
 
  901           evalGradDivu(du1, du2, du3,  x, y, z);
 
  902           rhsDatah(0,nPt,0) = du1;
 
  903           rhsDatah(0,nPt,1) = du2;
 
  904           rhsDatah(0,nPt,2) = du3;
 
  908       fst::integrate<double>(gC, rhsDatag, hexCurlsTransformedWeighted,
 
  912       fst::integrate<double>(hC, rhsDatah, hexCValsTransformedWeighted,
 
  916       fst::applyFieldSigns<double>(gC, hexEdgeSigns);
 
  917       fst::applyFieldSigns<double>(hC, hexEdgeSigns);
 
  921       for (
int i = 0; i < numFacesPerElem; i++){
 
  922         if (faceOnBoundary(elemToFace(k,i))){
 
  925             CellTools::mapToReferenceSubcell(refGaussPoints,
 
  930            hexHCurlBasis.
getValues(worksetCVals, refGaussPoints, OPERATOR_VALUE);
 
  933            CellTools::setJacobian(worksetJacobians,
 
  936            CellTools::setJacobianInv(worksetJacobInv, worksetJacobians );
 
  939             fst::HCURLtransformVALUE<double>(worksetCValsTransformed, worksetJacobInv,
 
  943             CellTools::mapToPhysicalFrame(worksetGaussPoints,
 
  948             CellTools::getPhysicalFaceNormals(worksetFaceN,
 
  953             for(
int nPt = 0; nPt < numFacePoints; nPt++){
 
  954                 for (
int dim = 0; dim < spaceDim; dim++){
 
  955                    worksetFaceNweighted(0,nPt,dim) = worksetFaceN(0,nPt,dim) * paramGaussWeights(nPt);
 
  959             fst::dotMultiplyDataField<double>(worksetFieldDotNormal, worksetFaceNweighted, 
 
  960                                                worksetCValsTransformed);
 
  963            for(
int nPt = 0; nPt < numFacePoints; nPt++){
 
  965              double x = worksetGaussPoints(0, nPt, 0);
 
  966              double y = worksetGaussPoints(0, nPt, 1);
 
  967              double z = worksetGaussPoints(0, nPt, 2);
 
  969              divuFace(0,nPt)=evalDivu(x, y, z);
 
  973           fst::integrate<double>(hCBoundary, divuFace, worksetFieldDotNormal,
 
  977            fst::applyFieldSigns<double>(hCBoundary, hexEdgeSigns);
 
  980             for (
int nF = 0; nF < numFieldsC; nF++){
 
  981                 hC(0,nF) = hC(0,nF) - hCBoundary(0,nF);
 
  989      for (
int row = 0; row < numFieldsC; row++){
 
  990            int rowIndex = elemToEdge(k,row);
 
  991            double val = gC(0,row)-hC(0,row);
 
  992            rhsC.SumIntoGlobalValues(1, &rowIndex, &val);
 
 1003    MassG.GlobalAssemble();  MassG.FillComplete();
 
 1004    MassC.GlobalAssemble();  MassC.FillComplete();
 
 1005    StiffC.GlobalAssemble(); StiffC.FillComplete();
 
 1006    rhsC.GlobalAssemble();
 
 1007    DGrad.GlobalAssemble(); DGrad.FillComplete(MassG.RowMap(),MassC.RowMap());
 
 1011    Epetra_CrsMatrix MassGinv(Copy,MassG.RowMap(),MassG.RowMap(),1);
 
 1012    Epetra_Vector DiagG(MassG.RowMap());
 
 1014    DiagG.PutScalar(1.0);
 
 1015    MassG.Multiply(
false,DiagG,DiagG);
 
 1016    for(
int i=0; i<DiagG.MyLength(); i++) {
 
 1017      DiagG[i]=1.0/DiagG[i];
 
 1019    for(
int i=0; i<DiagG.MyLength(); i++) {
 
 1020      int CID=MassG.GCID(i);
 
 1021      MassGinv.InsertGlobalValues(MassG.GRID(i),1,&(DiagG[i]),&CID);
 
 1023    MassGinv.FillComplete();
 
 1026    for(
int i=0;i<numNodes;i++) {
 
 1027      if (nodeOnBoundary(i)){
 
 1029       MassGinv.ReplaceGlobalValues(i,1,&val,&i);
 
 1038    for (
int row = 0; row<numEdges; row++){
 
 1039       MassC.ExtractMyRowView(row,numEntries,values,cols);
 
 1040         for (
int i=0; i<numEntries; i++){
 
 1041            if (edgeOnBoundary(cols[i])) {
 
 1045       StiffC.ExtractMyRowView(row,numEntries,values,cols);
 
 1046         for (
int i=0; i<numEntries; i++){
 
 1047            if (edgeOnBoundary(cols[i])) {
 
 1052    for (
int row = 0; row<numEdges; row++){
 
 1053        if (edgeOnBoundary(row)) {
 
 1055           StiffC.ExtractMyRowView(row,numEntries,values,cols);
 
 1056           for (
int i=0; i<numEntries; i++){
 
 1059           MassC.ExtractMyRowView(row,numEntries,values,cols);
 
 1060           for (
int i=0; i<numEntries; i++){
 
 1065          StiffC.ReplaceGlobalValues(1, &rowindex, 1, &rowindex, &val);
 
 1072    EpetraExt::RowMatrixToMatlabFile(
"mag_m0inv_matrix.dat",MassGinv);
 
 1073    EpetraExt::RowMatrixToMatlabFile(
"mag_m1_matrix.dat",MassC);
 
 1074    EpetraExt::RowMatrixToMatlabFile(
"mag_k1_matrix.dat",StiffC);
 
 1075    EpetraExt::RowMatrixToMatlabFile(
"mag_t0_matrix.dat",DGrad);
 
 1076    EpetraExt::MultiVectorToMatrixMarketFile(
"mag_rhs1_vector.dat",rhsC,0,0,
false);
 
 1081    std::cout << 
"End Result: TEST PASSED\n";
 
 1084  std::cout.copyfmt(oldFormatState);
 
 1090  int evalu(
double & uExact0, 
double & uExact1, 
double & uExact2, 
double & x, 
double & y, 
double & z)
 
 1093     uExact0 = exp(x+y+z)*(y+1.0)*(y-1.0)*(z+1.0)*(z-1.0);
 
 1094     uExact1 = exp(x+y+z)*(x+1.0)*(x-1.0)*(z+1.0)*(z-1.0);
 
 1095     uExact2 = exp(x+y+z)*(x+1.0)*(x-1.0)*(y+1.0)*(y-1.0);
 
 1118  double evalDivu(
double & x, 
double & y, 
double & z)
 
 1121     double divu = exp(x+y+z)*(y+1.0)*(y-1.0)*(z+1.0)*(z-1.0)
 
 1122                  + exp(x+y+z)*(x+1.0)*(x-1.0)*(z+1.0)*(z-1.0)
 
 1123                  + exp(x+y+z)*(x+1.0)*(x-1.0)*(y+1.0)*(y-1.0);
 
 1143  int evalCurlu(
double & curlu0, 
double & curlu1, 
double & curlu2, 
double & x, 
double & y, 
double & z)
 
 1147     double duxdy = exp(x+y+z)*(z*z-1.0)*(y*y+2.0*y-1.0);
 
 1148     double duxdz = exp(x+y+z)*(y*y-1.0)*(z*z+2.0*z-1.0);
 
 1149     double duydx = exp(x+y+z)*(z*z-1.0)*(x*x+2.0*x-1.0);
 
 1150     double duydz = exp(x+y+z)*(x*x-1.0)*(z*z+2.0*z-1.0);
 
 1151     double duzdx = exp(x+y+z)*(y*y-1.0)*(x*x+2.0*x-1.0);
 
 1152     double duzdy = exp(x+y+z)*(x*x-1.0)*(y*y+2.0*y-1.0);
 
 1180    curlu0 = duzdy - duydz;
 
 1181    curlu1 = duxdz - duzdx;
 
 1182    curlu2 = duydx - duxdy;
 
 1188  int evalGradDivu(
double & gradDivu0, 
double & gradDivu1, 
double & gradDivu2, 
double & x, 
double & y, 
double & z)
 
 1192     gradDivu0 = exp(x+y+z)*((y*y-1.0)*(z*z-1.0)+(x*x+2.0*x-1.0)*(z*z-1.0)+(x*x+2.0*x-1.0)*(y*y-1.0));
 
 1193     gradDivu1 = exp(x+y+z)*((y*y+2.0*y-1.0)*(z*z-1.0)+(x*x-1.0)*(z*z-1.0)+(x*x-1.0)*(y*y+2.0*y-1.0));
 
 1194     gradDivu2 = exp(x+y+z)*((y*y-1.0)*(z*z+2.0*z-1.0)+(x*x-1.0)*(z*z+2.0*z-1.0)+(x*x-1.0)*(y*y-1.0));
 
virtual int getCardinality() const 
Returns cardinality of the basis. 
Header file for utility class to provide multidimensional containers. 
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const 
Evaluation of a FEM basis on a reference Hexahedron cell. 
Header file for the abstract base class Intrepid::DefaultCubatureFactory. 
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Hexahedron cell...
Implementation of a templated lexicographical container for a multi-indexed scalar quantity...
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. 
Header file for the Intrepid::HCURL_HEX_I1_FEM class.