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