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