|
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 00080 // Intrepid includes 00081 #include "Intrepid_FunctionSpaceTools.hpp" 00082 #include "Intrepid_FieldContainer.hpp" 00083 #include "Intrepid_CellTools.hpp" 00084 #include "Intrepid_ArrayTools.hpp" 00085 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp" 00086 #include "Intrepid_RealSpaceTools.hpp" 00087 #include "Intrepid_DefaultCubatureFactory.hpp" 00088 #include "Intrepid_Utils.hpp" 00089 00090 // Epetra includes 00091 #include "Epetra_Time.h" 00092 #include "Epetra_Map.h" 00093 #include "Epetra_FECrsMatrix.h" 00094 #include "Epetra_FEVector.h" 00095 #include "Epetra_SerialComm.h" 00096 00097 // Teuchos includes 00098 #include "Teuchos_oblackholestream.hpp" 00099 #include "Teuchos_RCP.hpp" 00100 #include "Teuchos_BLAS.hpp" 00101 00102 // Shards includes 00103 #include "Shards_CellTopology.hpp" 00104 00105 // EpetraExt includes 00106 #include "EpetraExt_RowMatrixOut.h" 00107 #include "EpetraExt_MultiVectorOut.h" 00108 00109 using namespace std; 00110 using namespace Intrepid; 00111 00112 // Functions to evaluate exact solution and derivatives 00113 double evalu(double & x, double & y, double & z); 00114 int evalGradu(double & x, double & y, double & z, double & gradu1, double & gradu2, double & gradu3); 00115 double evalDivGradu(double & x, double & y, double & z); 00116 00117 int main(int argc, char *argv[]) { 00118 00119 //Check number of arguments 00120 if (argc < 4) { 00121 std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n"; 00122 std::cout <<"Usage:\n\n"; 00123 std::cout <<" ./Intrepid_example_Drivers_Example_03.exe NX NY NZ verbose\n\n"; 00124 std::cout <<" where \n"; 00125 std::cout <<" int NX - num intervals in x direction (assumed box domain, 0,1) \n"; 00126 std::cout <<" int NY - num intervals in y direction (assumed box domain, 0,1) \n"; 00127 std::cout <<" int NZ - num intervals in z direction (assumed box domain, 0,1) \n"; 00128 std::cout <<" verbose (optional) - any character, indicates verbose output \n\n"; 00129 exit(1); 00130 } 00131 00132 // This little trick lets us print to std::cout only if 00133 // a (dummy) command-line argument is provided. 00134 int iprint = argc - 1; 00135 Teuchos::RCP<std::ostream> outStream; 00136 Teuchos::oblackholestream bhs; // outputs nothing 00137 if (iprint > 3) 00138 outStream = Teuchos::rcp(&std::cout, false); 00139 else 00140 outStream = Teuchos::rcp(&bhs, false); 00141 00142 // Save the format state of the original std::cout. 00143 Teuchos::oblackholestream oldFormatState; 00144 oldFormatState.copyfmt(std::cout); 00145 00146 *outStream \ 00147 << "===============================================================================\n" \ 00148 << "| |\n" \ 00149 << "| Example: Generate Stiffness Matrix and Right Hand Side Vector for |\n" \ 00150 << "| Poisson Equation on Hexahedral Mesh |\n" \ 00151 << "| |\n" \ 00152 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00153 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00154 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00155 << "| |\n" \ 00156 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00157 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00158 << "| |\n" \ 00159 << "===============================================================================\n"; 00160 00161 00162 // ************************************ GET INPUTS ************************************** 00163 00164 int NX = atoi(argv[1]); // num intervals in x direction (assumed box domain, 0,1) 00165 int NY = atoi(argv[2]); // num intervals in y direction (assumed box domain, 0,1) 00166 int NZ = atoi(argv[3]); // num intervals in z direction (assumed box domain, 0,1) 00167 00168 // *********************************** CELL TOPOLOGY ********************************** 00169 00170 // Get cell topology for base hexahedron 00171 typedef shards::CellTopology CellTopology; 00172 CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() ); 00173 00174 // Get dimensions 00175 int numNodesPerElem = hex_8.getNodeCount(); 00176 int spaceDim = hex_8.getDimension(); 00177 00178 // *********************************** GENERATE MESH ************************************ 00179 00180 *outStream << "Generating mesh ... \n\n"; 00181 00182 *outStream << " NX" << " NY" << " NZ\n"; 00183 *outStream << std::setw(5) << NX << 00184 std::setw(5) << NY << 00185 std::setw(5) << NZ << "\n\n"; 00186 00187 // Print mesh information 00188 int numElems = NX*NY*NZ; 00189 int numNodes = (NX+1)*(NY+1)*(NZ+1); 00190 *outStream << " Number of Elements: " << numElems << " \n"; 00191 *outStream << " Number of Nodes: " << numNodes << " \n\n"; 00192 00193 // Cube 00194 double leftX = 0.0, rightX = 1.0; 00195 double leftY = 0.0, rightY = 1.0; 00196 double leftZ = 0.0, rightZ = 1.0; 00197 00198 // Mesh spacing 00199 double hx = (rightX-leftX)/((double)NX); 00200 double hy = (rightY-leftY)/((double)NY); 00201 double hz = (rightZ-leftZ)/((double)NZ); 00202 00203 // Get nodal coordinates 00204 FieldContainer<double> nodeCoord(numNodes, spaceDim); 00205 FieldContainer<int> nodeOnBoundary(numNodes); 00206 int inode = 0; 00207 for (int k=0; k<NZ+1; k++) { 00208 for (int j=0; j<NY+1; j++) { 00209 for (int i=0; i<NX+1; i++) { 00210 nodeCoord(inode,0) = leftX + (double)i*hx; 00211 nodeCoord(inode,1) = leftY + (double)j*hy; 00212 nodeCoord(inode,2) = leftZ + (double)k*hz; 00213 if (k==0 || j==0 || i==0 || k==NZ || j==NY || i==NX){ 00214 nodeOnBoundary(inode)=1; 00215 } 00216 else { 00217 nodeOnBoundary(inode)=0; 00218 } 00219 inode++; 00220 } 00221 } 00222 } 00223 #define DUMP_DATA 00224 #ifdef DUMP_DATA 00225 // Print nodal coords 00226 ofstream fcoordout("coords.dat"); 00227 for (int i=0; i<numNodes; i++) { 00228 fcoordout << nodeCoord(i,0) <<" "; 00229 fcoordout << nodeCoord(i,1) <<" "; 00230 fcoordout << nodeCoord(i,2) <<"\n"; 00231 } 00232 fcoordout.close(); 00233 #endif 00234 00235 00236 // Element to Node map 00237 FieldContainer<int> elemToNode(numElems, numNodesPerElem); 00238 int ielem = 0; 00239 for (int k=0; k<NZ; k++) { 00240 for (int j=0; j<NY; j++) { 00241 for (int i=0; i<NX; i++) { 00242 elemToNode(ielem,0) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i; 00243 elemToNode(ielem,1) = (NY + 1)*(NX + 1)*k + (NX + 1)*j + i + 1; 00244 elemToNode(ielem,2) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i + 1; 00245 elemToNode(ielem,3) = (NY + 1)*(NX + 1)*k + (NX + 1)*(j + 1) + i; 00246 elemToNode(ielem,4) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i; 00247 elemToNode(ielem,5) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*j + i + 1; 00248 elemToNode(ielem,6) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i + 1; 00249 elemToNode(ielem,7) = (NY + 1)*(NX + 1)*(k + 1) + (NX + 1)*(j + 1) + i; 00250 ielem++; 00251 } 00252 } 00253 } 00254 #ifdef DUMP_DATA 00255 // Output connectivity 00256 ofstream fe2nout("elem2node.dat"); 00257 for (int k=0; k<NZ; k++) { 00258 for (int j=0; j<NY; j++) { 00259 for (int i=0; i<NX; i++) { 00260 int ielem = i + j * NX + k * NX * NY; 00261 for (int m=0; m<numNodesPerElem; m++){ 00262 fe2nout << elemToNode(ielem,m) <<" "; 00263 } 00264 fe2nout <<"\n"; 00265 } 00266 } 00267 } 00268 fe2nout.close(); 00269 #endif 00270 00271 00272 // ************************************ CUBATURE ************************************** 00273 00274 *outStream << "Getting cubature ... \n\n"; 00275 00276 // Get numerical integration points and weights 00277 DefaultCubatureFactory<double> cubFactory; 00278 int cubDegree = 2; 00279 Teuchos::RCP<Cubature<double> > hexCub = cubFactory.create(hex_8, cubDegree); 00280 00281 int cubDim = hexCub->getDimension(); 00282 int numCubPoints = hexCub->getNumPoints(); 00283 00284 FieldContainer<double> cubPoints(numCubPoints, cubDim); 00285 FieldContainer<double> cubWeights(numCubPoints); 00286 00287 hexCub->getCubature(cubPoints, cubWeights); 00288 00289 00290 // ************************************** BASIS *************************************** 00291 00292 *outStream << "Getting basis ... \n\n"; 00293 00294 // Define basis 00295 Basis_HGRAD_HEX_C1_FEM<double, FieldContainer<double> > hexHGradBasis; 00296 int numFieldsG = hexHGradBasis.getCardinality(); 00297 FieldContainer<double> hexGVals(numFieldsG, numCubPoints); 00298 FieldContainer<double> hexGrads(numFieldsG, numCubPoints, spaceDim); 00299 00300 // Evaluate basis values and gradients at cubature points 00301 hexHGradBasis.getValues(hexGVals, cubPoints, OPERATOR_VALUE); 00302 hexHGradBasis.getValues(hexGrads, cubPoints, OPERATOR_GRAD); 00303 00304 00305 // ******** LOOP OVER ELEMENTS TO CREATE LOCAL STIFFNESS MATRIX ************* 00306 00307 *outStream << "Building stiffness matrix and right hand side ... \n\n"; 00308 00309 // Settings and data structures for mass and stiffness matrices 00310 typedef CellTools<double> CellTools; 00311 typedef FunctionSpaceTools fst; 00312 int numCells = 1; 00313 00314 // Container for nodes 00315 FieldContainer<double> hexNodes(numCells, numNodesPerElem, spaceDim); 00316 // Containers for Jacobian 00317 FieldContainer<double> hexJacobian(numCells, numCubPoints, spaceDim, spaceDim); 00318 FieldContainer<double> hexJacobInv(numCells, numCubPoints, spaceDim, spaceDim); 00319 FieldContainer<double> hexJacobDet(numCells, numCubPoints); 00320 // Containers for element HGRAD stiffness matrix 00321 FieldContainer<double> localStiffMatrix(numCells, numFieldsG, numFieldsG); 00322 FieldContainer<double> weightedMeasure(numCells, numCubPoints); 00323 FieldContainer<double> hexGradsTransformed(numCells, numFieldsG, numCubPoints, spaceDim); 00324 FieldContainer<double> hexGradsTransformedWeighted(numCells, numFieldsG, numCubPoints, spaceDim); 00325 // Containers for right hand side vectors 00326 FieldContainer<double> rhsData(numCells, numCubPoints); 00327 FieldContainer<double> localRHS(numCells, numFieldsG); 00328 FieldContainer<double> hexGValsTransformed(numCells, numFieldsG, numCubPoints); 00329 FieldContainer<double> hexGValsTransformedWeighted(numCells, numFieldsG, numCubPoints); 00330 // Container for cubature points in physical space 00331 FieldContainer<double> physCubPoints(numCells, numCubPoints, cubDim); 00332 00333 // Global arrays in Epetra format 00334 Epetra_SerialComm Comm; 00335 Epetra_Map globalMapG(numNodes, 0, Comm); 00336 Epetra_FECrsMatrix StiffMatrix(Copy, globalMapG, numFieldsG); 00337 Epetra_FEVector rhs(globalMapG); 00338 00339 // *** Element loop *** 00340 for (int k=0; k<numElems; k++) { 00341 00342 // Physical cell coordinates 00343 for (int i=0; i<numNodesPerElem; i++) { 00344 hexNodes(0,i,0) = nodeCoord(elemToNode(k,i),0); 00345 hexNodes(0,i,1) = nodeCoord(elemToNode(k,i),1); 00346 hexNodes(0,i,2) = nodeCoord(elemToNode(k,i),2); 00347 } 00348 00349 // Compute cell Jacobians, their inverses and their determinants 00350 CellTools::setJacobian(hexJacobian, cubPoints, hexNodes, hex_8); 00351 CellTools::setJacobianInv(hexJacobInv, hexJacobian ); 00352 CellTools::setJacobianDet(hexJacobDet, hexJacobian ); 00353 00354 // ************************** Compute element HGrad stiffness matrices ******************************* 00355 00356 // transform to physical coordinates 00357 fst::HGRADtransformGRAD<double>(hexGradsTransformed, hexJacobInv, hexGrads); 00358 00359 // compute weighted measure 00360 fst::computeCellMeasure<double>(weightedMeasure, hexJacobDet, cubWeights); 00361 00362 // multiply values with weighted measure 00363 fst::multiplyMeasure<double>(hexGradsTransformedWeighted, 00364 weightedMeasure, hexGradsTransformed); 00365 00366 // integrate to compute element stiffness matrix 00367 fst::integrate<double>(localStiffMatrix, 00368 hexGradsTransformed, hexGradsTransformedWeighted, COMP_BLAS); 00369 00370 // assemble into global matrix 00371 for (int row = 0; row < numFieldsG; row++){ 00372 for (int col = 0; col < numFieldsG; col++){ 00373 int rowIndex = elemToNode(k,row); 00374 int colIndex = elemToNode(k,col); 00375 double val = localStiffMatrix(0,row,col); 00376 StiffMatrix.InsertGlobalValues(1, &rowIndex, 1, &colIndex, &val); 00377 } 00378 } 00379 00380 // ******************************* Build right hand side ************************************ 00381 00382 // transform integration points to physical points 00383 CellTools::mapToPhysicalFrame(physCubPoints, cubPoints, hexNodes, hex_8); 00384 00385 // evaluate right hand side function at physical points 00386 for (int nPt = 0; nPt < numCubPoints; nPt++){ 00387 00388 double x = physCubPoints(0,nPt,0); 00389 double y = physCubPoints(0,nPt,1); 00390 double z = physCubPoints(0,nPt,2); 00391 00392 rhsData(0,nPt) = evalDivGradu(x, y, z); 00393 } 00394 00395 // transform basis values to physical coordinates 00396 fst::HGRADtransformVALUE<double>(hexGValsTransformed, hexGVals); 00397 00398 // multiply values with weighted measure 00399 fst::multiplyMeasure<double>(hexGValsTransformedWeighted, 00400 weightedMeasure, hexGValsTransformed); 00401 00402 // integrate rhs term 00403 fst::integrate<double>(localRHS, rhsData, hexGValsTransformedWeighted, 00404 COMP_BLAS); 00405 00406 // assemble into global vector 00407 for (int row = 0; row < numFieldsG; row++){ 00408 int rowIndex = elemToNode(k,row); 00409 double val = -localRHS(0,row); 00410 rhs.SumIntoGlobalValues(1, &rowIndex, &val); 00411 } 00412 00413 } // *** end element loop *** 00414 00415 00416 // Assemble global matrices 00417 StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete(); 00418 rhs.GlobalAssemble(); 00419 00420 00421 // Adjust stiffness matrix and rhs based on boundary conditions 00422 for (int row = 0; row<numNodes; row++){ 00423 if (nodeOnBoundary(row)) { 00424 int rowindex = row; 00425 for (int col=0; col<numNodes; col++){ 00426 double val = 0.0; 00427 int colindex = col; 00428 StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &colindex, &val); 00429 } 00430 double val = 1.0; 00431 StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &rowindex, &val); 00432 val = 0.0; 00433 rhs.ReplaceGlobalValues(1, &rowindex, &val); 00434 } 00435 } 00436 00437 #ifdef DUMP_DATA 00438 // Dump matrices to disk 00439 EpetraExt::RowMatrixToMatlabFile("stiff_matrix.dat",StiffMatrix); 00440 EpetraExt::MultiVectorToMatrixMarketFile("rhs_vector.dat",rhs,0,0,false); 00441 #endif 00442 00443 std::cout << "End Result: TEST PASSED\n"; 00444 00445 // reset format state of std::cout 00446 std::cout.copyfmt(oldFormatState); 00447 00448 return 0; 00449 } 00450 00451 00452 // Calculates value of exact solution u 00453 double evalu(double & x, double & y, double & z) 00454 { 00455 /* 00456 // function1 00457 double exactu = sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z); 00458 */ 00459 00460 // function2 00461 double exactu = sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z); 00462 00463 return exactu; 00464 } 00465 00466 // Calculates gradient of exact solution u 00467 int evalGradu(double & x, double & y, double & z, double & gradu1, double & gradu2, double & gradu3) 00468 { 00469 /* 00470 // function 1 00471 gradu1 = M_PI*cos(M_PI*x)*sin(M_PI*y)*sin(M_PI*z); 00472 gradu2 = M_PI*sin(M_PI*x)*cos(M_PI*y)*sin(M_PI*z); 00473 gradu3 = M_PI*sin(M_PI*x)*sin(M_PI*y)*cos(M_PI*z); 00474 */ 00475 00476 // function2 00477 gradu1 = (M_PI*cos(M_PI*x)+sin(M_PI*x)) 00478 *sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z); 00479 gradu2 = (M_PI*cos(M_PI*y)+sin(M_PI*y)) 00480 *sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z); 00481 gradu3 = (M_PI*cos(M_PI*z)+sin(M_PI*z)) 00482 *sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z); 00483 00484 return 0; 00485 } 00486 00487 // Calculates Laplacian of exact solution u 00488 double evalDivGradu(double & x, double & y, double & z) 00489 { 00490 /* 00491 // function 1 00492 double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z); 00493 */ 00494 00495 // function 2 00496 double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z) 00497 + 2.0*M_PI*cos(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z) 00498 + 2.0*M_PI*cos(M_PI*y)*sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z) 00499 + 2.0*M_PI*cos(M_PI*z)*sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z) 00500 + 3.0*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z); 00501 00502 return divGradu; 00503 } 00504
1.7.6.1