|
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 00069 // Intrepid includes 00070 #include "Intrepid_FunctionSpaceTools.hpp" 00071 #include "Intrepid_FieldContainer.hpp" 00072 #include "Intrepid_CellTools.hpp" 00073 #include "Intrepid_CubaturePolylib.hpp" 00074 //#include "Intrepid_ArrayTools.hpp" 00075 #include "Intrepid_HGRAD_LINE_Cn_FEM.hpp" 00076 //#include "Intrepid_RealSpaceTools.hpp" 00077 #include "Intrepid_Utils.hpp" 00078 00079 // Epetra includes 00080 #include "Epetra_Time.h" 00081 #include "Epetra_Map.h" 00082 #include "Epetra_FEVector.h" 00083 #include "Epetra_SerialComm.h" 00084 00085 // Teuchos includes 00086 #include "Teuchos_oblackholestream.hpp" 00087 #include "Teuchos_RCP.hpp" 00088 //#include "Teuchos_BLAS.hpp" 00089 //#include "Teuchos_BLAS_types.hpp" 00090 00091 // Shards includes 00092 #include "Shards_CellTopology.hpp" 00093 00094 // EpetraExt includes 00095 #include "EpetraExt_MultiVectorOut.h" 00096 00097 using namespace std; 00098 using namespace Intrepid; 00099 00100 int main(int argc, char *argv[]) { 00101 00102 //Check number of arguments 00103 if (argc < 4) { 00104 std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n"; 00105 std::cout <<"Usage:\n\n"; 00106 std::cout <<" ./Intrepid_example_Drivers_Example_14.exe deg NX NY NZ verbose\n\n"; 00107 std::cout <<" where \n"; 00108 std::cout <<" int deg - polynomial degree to be used (assumed >= 1) \n"; 00109 std::cout <<" int NX - num intervals in x direction (assumed box domain, 0,1) \n"; 00110 std::cout <<" int NY - num intervals in y direction (assumed box domain, 0,1) \n"; 00111 std::cout <<" int NZ - num intervals in y direction (assumed box domain, 0,1) \n"; 00112 std::cout <<" verbose (optional) - any character, indicates verbose output \n\n"; 00113 exit(1); 00114 } 00115 00116 // This little trick lets us print to std::cout only if 00117 // a (dummy) command-line argument is provided. 00118 int iprint = argc - 1; 00119 Teuchos::RCP<std::ostream> outStream; 00120 Teuchos::oblackholestream bhs; // outputs nothing 00121 if (iprint > 2) 00122 outStream = Teuchos::rcp(&std::cout, false); 00123 else 00124 outStream = Teuchos::rcp(&bhs, false); 00125 00126 // Save the format state of the original std::cout. 00127 Teuchos::oblackholestream oldFormatState; 00128 oldFormatState.copyfmt(std::cout); 00129 00130 *outStream \ 00131 << "===============================================================================\n" \ 00132 << "| |\n" \ 00133 << "| Example: Apply Stiffness Matrix for |\n" \ 00134 << "| Poisson Equation on Hexahedral Mesh |\n" \ 00135 << "| |\n" \ 00136 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00137 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00138 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00139 << "| |\n" \ 00140 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00141 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00142 << "| |\n" \ 00143 << "===============================================================================\n"; 00144 00145 00146 // ************************************ GET INPUTS ************************************** 00147 00148 int deg = atoi(argv[1]); // polynomial degree to use 00149 int NX = atoi(argv[2]); // num intervals in x direction (assumed box domain, 0,1) 00150 int NY = atoi(argv[3]); // num intervals in y direction (assumed box domain, 0,1) 00151 int NZ = atoi(argv[4]); // num intervals in y direction (assumed box domain, 0,1) 00152 00153 00154 // *********************************** CELL TOPOLOGY ********************************** 00155 00156 // Get cell topology for base hexahedron 00157 typedef shards::CellTopology CellTopology; 00158 CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() ); 00159 00160 // Get dimensions 00161 int numNodesPerElem = hex_8.getNodeCount(); 00162 int spaceDim = hex_8.getDimension(); 00163 00164 // *********************************** GENERATE MESH ************************************ 00165 00166 *outStream << "Generating mesh ... \n\n"; 00167 00168 *outStream << " NX" << " NY" << " NZ\n"; 00169 *outStream << std::setw(5) << NX << 00170 std::setw(5) << NY << std::setw(5) << NZ << "\n\n"; 00171 00172 // Print mesh information 00173 int numElems = NX*NY*NZ; 00174 int numNodes = (NX+1)*(NY+1)*(NZ+1); 00175 *outStream << " Number of Elements: " << numElems << " \n"; 00176 *outStream << " Number of Nodes: " << numNodes << " \n\n"; 00177 00178 // Cube 00179 double leftX = 0.0, rightX = 1.0; 00180 double leftY = 0.0, rightY = 1.0; 00181 double leftZ = 0.0, rightZ = 1.0; 00182 00183 // Mesh spacing 00184 double hx = (rightX-leftX)/((double)NX); 00185 double hy = (rightY-leftY)/((double)NY); 00186 double hz = (rightZ-leftZ)/((double)NZ); 00187 00188 // Get nodal coordinates 00189 FieldContainer<double> nodeCoord(numNodes, spaceDim); 00190 FieldContainer<int> nodeOnBoundary(numNodes); 00191 int inode = 0; 00192 for (int k=0; k<NZ+1; k++) 00193 { 00194 for (int j=0; j<NY+1; j++) 00195 { 00196 for (int i=0; i<NX+1; i++) 00197 { 00198 nodeCoord(inode,0) = leftX + (double)i*hx; 00199 nodeCoord(inode,1) = leftY + (double)j*hy; 00200 nodeCoord(inode,2) = leftZ + (double)k*hz; 00201 if (k==0 || k==NZ || j==0 || i==0 || j==NY || i==NX) 00202 { 00203 nodeOnBoundary(inode)=1; 00204 } 00205 else 00206 { 00207 nodeOnBoundary(inode)=0; 00208 } 00209 inode++; 00210 } 00211 } 00212 } 00213 #define DUMP_DATA 00214 #ifdef DUMP_DATA 00215 // Print nodal coords 00216 ofstream fcoordout("coords.dat"); 00217 for (int i=0; i<numNodes; i++) { 00218 fcoordout << nodeCoord(i,0) <<" "; 00219 fcoordout << nodeCoord(i,1) <<" "; 00220 fcoordout << nodeCoord(i,2) <<"\n"; 00221 } 00222 fcoordout.close(); 00223 #endif 00224 00225 00226 // Element to Node map 00227 // We'll keep it around, but this is only the DOFMap if you are in the lowest order case. 00228 FieldContainer<int> elemToNode(numElems, numNodesPerElem); 00229 int ielem = 0; 00230 for (int k=0; k<NZ; k++) 00231 { 00232 for (int j=0; j<NY; j++) 00233 { 00234 for (int i=0; i<NX; i++) 00235 { 00236 elemToNode(ielem,0) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i; 00237 elemToNode(ielem,1) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1; 00238 elemToNode(ielem,2) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1; 00239 elemToNode(ielem,3) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i; 00240 elemToNode(ielem,4) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i; 00241 elemToNode(ielem,5) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1; 00242 elemToNode(ielem,6) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1; 00243 elemToNode(ielem,7) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i; 00244 ielem++; 00245 } 00246 } 00247 } 00248 #ifdef DUMP_DATA 00249 // Output connectivity 00250 ofstream fe2nout("elem2node.dat"); 00251 for (int k=0;k<NZ;k++) 00252 { 00253 for (int j=0; j<NY; j++) 00254 { 00255 for (int i=0; i<NX; i++) 00256 { 00257 int ielem = i + j * NX + k * NY * NY; 00258 for (int m=0; m<numNodesPerElem; m++) 00259 { 00260 fe2nout << elemToNode(ielem,m) <<" "; 00261 } 00262 fe2nout <<"\n"; 00263 } 00264 } 00265 } 00266 fe2nout.close(); 00267 #endif 00268 00269 // ********************************* 1-D CUBATURE AND BASIS *********************************** 00270 *outStream << "Getting cubature and basis ... \n\n"; 00271 00272 // Get numerical integration points and weights 00273 // I only need this on the line since I'm doing tensor products 00274 Teuchos::RCP<Cubature<double,FieldContainer<double>,FieldContainer<double> > > glcub 00275 = Teuchos::rcp(new CubaturePolylib<double,FieldContainer<double>,FieldContainer<double> >(2*deg-1,PL_GAUSS_LOBATTO) ); 00276 00277 const int numCubPoints = glcub->getNumPoints(); 00278 00279 FieldContainer<double> cubPoints1D(numCubPoints, 1); 00280 FieldContainer<double> cubWeights1D(numCubPoints); 00281 00282 glcub->getCubature(cubPoints1D,cubWeights1D); 00283 // Define basis: I only need this on the line also 00284 Basis_HGRAD_LINE_Cn_FEM<double, FieldContainer<double> > lineHGradBasis(deg,POINTTYPE_SPECTRAL); 00285 int numLineFieldsG = lineHGradBasis.getCardinality(); 00286 FieldContainer<double> lineGrads(numLineFieldsG, numCubPoints, 1); 00287 00288 // Evaluate basis values and gradients at cubature points 00289 lineHGradBasis.getValues(lineGrads, cubPoints1D, OPERATOR_GRAD); 00290 00291 00292 // ********************************* 3-D LOCAL-TO-GLOBAL MAPPING ******************************* 00293 FieldContainer<int> ltgMapping(numElems,numLineFieldsG*numLineFieldsG*numLineFieldsG); 00294 const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1); 00295 ielem=0; 00296 for (int k=0;k<NZ;k++) 00297 { 00298 for (int j=0;j<NY;j++) 00299 { 00300 for (int i=0;i<NX;i++) 00301 { 00302 const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg; 00303 // loop over local dof on this cell 00304 int local_dof_cur=0; 00305 for (int kloc=0;kloc<=deg;kloc++) 00306 { 00307 for (int jloc=0;jloc<=deg;jloc++) 00308 { 00309 for (int iloc=0;iloc<=deg;iloc++) 00310 { 00311 ltgMapping(ielem,local_dof_cur) = start 00312 + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 ) 00313 + jloc * ( NX * deg + 1 ) 00314 + iloc; 00315 local_dof_cur++; 00316 } 00317 } 00318 } 00319 ielem++; 00320 } 00321 } 00322 } 00323 #ifdef DUMP_DATA 00324 // Output ltg mapping 00325 ielem = 0; 00326 ofstream ltgout("ltg.dat"); 00327 for (int k=0;k<NZ;k++) 00328 { 00329 for (int j=0; j<NY; j++) 00330 { 00331 for (int i=0; i<NX; i++) 00332 { 00333 int ielem = i + j * NX + k * NX * NY; 00334 for (int m=0; m<numLineFieldsG*numLineFieldsG*numLineFieldsG; m++) 00335 { 00336 ltgout << ltgMapping(ielem,m) <<" "; 00337 } 00338 ltgout <<"\n"; 00339 } 00340 } 00341 } 00342 ltgout.close(); 00343 #endif 00344 00345 // ********** DECLARE GLOBAL OBJECTS ************* 00346 Epetra_SerialComm Comm; 00347 Epetra_Map globalMapG(numDOF, 0, Comm); 00348 Epetra_FEVector u(globalMapG); u.Random(); 00349 Epetra_FEVector Ku(globalMapG); 00350 00351 00352 // ************* MATRIX-FREE APPLICATION 00353 FieldContainer<double> uScattered(numElems,numLineFieldsG*numLineFieldsG*numLineFieldsG); 00354 FieldContainer<double> KuScattered(numElems,numLineFieldsG*numLineFieldsG*numLineFieldsG); 00355 00356 u.GlobalAssemble(); 00357 00358 Epetra_Time multTimer(Comm); 00359 Teuchos::BLAS<int,double> blas; 00360 Ku.PutScalar(0.0); 00361 Ku.GlobalAssemble(); 00362 00363 double *uVals = u[0]; 00364 double *KuVals = Ku[0]; 00365 00366 Epetra_Time scatterTimer(Comm); 00367 std::cout << "Scattering\n"; 00368 // Scatter 00369 for (int k=0; k<numElems; k++) 00370 { 00371 for (int i=0;i<numLineFieldsG*numLineFieldsG*numLineFieldsG;i++) 00372 { 00373 uScattered(k,i) = uVals[ltgMapping(k,i)]; 00374 } 00375 } 00376 00377 00378 const double scatterTime = scatterTimer.ElapsedTime(); 00379 00380 00381 00382 FieldContainer<double> Du(numLineFieldsG,numLineFieldsG,numLineFieldsG); 00383 00384 Epetra_Time localAppTimer(Comm); 00385 00386 uScattered.resize(numElems,numLineFieldsG,numLineFieldsG,numLineFieldsG); 00387 00388 00389 int cur; 00390 double hcur; 00391 00392 for (ielem=0;ielem<numElems;ielem++) 00393 { 00394 // X-COMPONENT OF ELEMENT LAPLACIAN 00395 00396 // zero out derivative 00397 for (int k=0;k<numLineFieldsG;k++) 00398 { 00399 for (int j=0;j<numLineFieldsG;j++) 00400 { 00401 for (int i=0;i<numLineFieldsG;i++) 00402 { 00403 Du(k,j,i) = 0.0; 00404 } 00405 } 00406 } 00407 00408 00409 // compute x derivative 00410 // this could be a very simple dgemm call 00411 for (int k=0;k<numLineFieldsG;k++) 00412 { 00413 for (int j=0;j<numLineFieldsG;j++) 00414 { 00415 for (int i=0;i<numLineFieldsG;i++) 00416 { 00417 for (int q=0;q<numLineFieldsG;q++) 00418 { 00419 Du(k,j,i) += uScattered(ielem,k,j,i) * lineGrads(i,q,0); 00420 } 00421 } 00422 } 00423 } 00424 00425 // integration loop for x derivative term 00426 cur = 0; 00427 hcur = hy * hz / hx; 00428 for (int k=0;k<numLineFieldsG;k++) 00429 { 00430 const double wt1 = hcur * cubWeights1D(k); 00431 for (int j=0;j<numLineFieldsG;j++) 00432 { 00433 const double wtstuff = wt1 * cubWeights1D(j); 00434 for (int i=0;i<numLineFieldsG;i++) 00435 { 00436 for (int q=0;q<numLineFieldsG;q++) 00437 { 00438 KuScattered(ielem,cur) += wtstuff 00439 * cubWeights1D(q) * Du(k,j,q) * lineGrads(i,q,0); 00440 } 00441 cur++; 00442 } 00443 } 00444 } 00445 00446 00447 // Y-COMPONENT OF ELEMENT LAPLACIAN 00448 00449 // zero out derivative 00450 for (int k=0;k<numLineFieldsG;k++) 00451 { 00452 for (int j=0;j<numLineFieldsG;j++) 00453 { 00454 for (int i=0;i<numLineFieldsG;i++) 00455 { 00456 Du(k,j,i) = 0.0; 00457 } 00458 } 00459 } 00460 00461 // compute y derivative 00462 for (int k=0;k<numLineFieldsG;k++) 00463 { 00464 for (int j=0;j<numLineFieldsG;j++) 00465 { 00466 for (int i=0;i<numLineFieldsG;i++) 00467 { 00468 for (int q=0;q<numLineFieldsG;q++) 00469 { 00470 Du(k,j,i) += uScattered(ielem,k,j,i) * lineGrads(j,q,0); 00471 } 00472 } 00473 } 00474 } 00475 00476 // integration loop for y derivative term 00477 cur = 0; 00478 hcur = hx * hz / hy; 00479 for (int k=0;k<numLineFieldsG;k++) 00480 { 00481 const double wt1 = hcur * cubWeights1D(k); 00482 for (int j=0;j<numLineFieldsG;j++) 00483 { 00484 for (int i=0;i<numLineFieldsG;i++) 00485 { 00486 const double wtstuff = cubWeights1D(i) * wt1; 00487 for (int q=0;q<numLineFieldsG;q++) 00488 { 00489 KuScattered(ielem,cur) += wtstuff * cubWeights1D(q) * Du(k,q,i) * lineGrads(j,q,0); 00490 } 00491 cur++; 00492 } 00493 } 00494 } 00495 00496 00497 // Z-COMPONENT OF ELEMENT LAPLACIAN 00498 00499 // zero out derivative 00500 for (int k=0;k<numLineFieldsG;k++) 00501 { 00502 for (int j=0;j<numLineFieldsG;j++) 00503 { 00504 for (int i=0;i<numLineFieldsG;i++) 00505 { 00506 Du(k,j,i) = 0.0; 00507 } 00508 } 00509 } 00510 00511 // compute z derivative 00512 for (int k=0;k<numLineFieldsG;k++) 00513 { 00514 for (int j=0;j<numLineFieldsG;j++) 00515 { 00516 for (int i=0;i<numLineFieldsG;i++) 00517 { 00518 for (int q=0;q<numLineFieldsG;q++) 00519 { 00520 Du(k,j,i) += uScattered(ielem,k,j,i) * lineGrads(k,q,0); 00521 } 00522 } 00523 } 00524 } 00525 00526 // integration loop for z derivative term. 00527 cur = 0; 00528 hcur = hx * hy / hz; 00529 for (int k=0;k<numLineFieldsG;k++) 00530 { 00531 for (int j=0;j<numLineFieldsG;j++) 00532 { 00533 const double wt1 = hcur * cubWeights1D(j); 00534 for (int i=0;i<numLineFieldsG;i++) 00535 { 00536 const double wtstuff = cubWeights1D(i) * wt1; 00537 for (int q=0;q<numLineFieldsG;q++) 00538 { 00539 KuScattered(ielem,cur) += wtstuff * cubWeights1D(q) * Du(q,j,i) * lineGrads(k,q,0); 00540 } 00541 cur++; 00542 } 00543 } 00544 } 00545 00546 } 00547 00548 const double localAppTime = localAppTimer.ElapsedTime(); 00549 00550 Epetra_Time gatherTimer(Comm); 00551 // Gather 00552 for (int k=0;k<numElems;k++) 00553 { 00554 for (int i=0;i<numLineFieldsG*numLineFieldsG*numLineFieldsG;i++) 00555 { 00556 KuVals[ltgMapping(k,i)] += KuScattered(k,i); 00557 } 00558 } 00559 00560 const double gatherTime = gatherTimer.ElapsedTime(); 00561 00562 00563 *outStream << "Time to scatter " << scatterTime << "\n"; 00564 *outStream << "Time for local application " << localAppTime << "\n"; 00565 *outStream << "Time to gather " << gatherTime << "\n"; 00566 *outStream << "Total matrix-free time " << scatterTime + localAppTime + gatherTime << "\n"; 00567 00568 00569 *outStream << "End Result: TEST PASSED\n"; 00570 00571 // reset format state of std::cout 00572 std::cout.copyfmt(oldFormatState); 00573 00574 return 0; 00575 } 00576
1.7.4