|
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 00068 // Intrepid includes 00069 #include "Intrepid_CubaturePolylib.hpp" 00070 #include "Intrepid_FunctionSpaceTools.hpp" 00071 #include "Intrepid_FieldContainer.hpp" 00072 #include "Intrepid_CellTools.hpp" 00073 #include "Intrepid_ArrayTools.hpp" 00074 #include "Intrepid_HGRAD_QUAD_Cn_FEM.hpp" 00075 #include "Intrepid_RealSpaceTools.hpp" 00076 #include "Intrepid_DefaultCubatureFactory.hpp" 00077 #include "Intrepid_Utils.hpp" 00078 00079 // Epetra includes 00080 #include "Epetra_Time.h" 00081 #include "Epetra_Map.h" 00082 #include "Epetra_FECrsMatrix.h" 00083 #include "Epetra_FEVector.h" 00084 #include "Epetra_SerialComm.h" 00085 00086 // Teuchos includes 00087 #include "Teuchos_oblackholestream.hpp" 00088 #include "Teuchos_RCP.hpp" 00089 #include "Teuchos_BLAS.hpp" 00090 00091 // Shards includes 00092 #include "Shards_CellTopology.hpp" 00093 00094 // EpetraExt includes 00095 #include "EpetraExt_RowMatrixOut.h" 00096 #include "EpetraExt_MultiVectorOut.h" 00097 00098 using namespace std; 00099 using namespace Intrepid; 00100 00101 int main(int argc, char *argv[]) { 00102 00103 //Check number of arguments 00104 if (argc < 4) { 00105 std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n"; 00106 std::cout <<"Usage:\n\n"; 00107 std::cout <<" ./Intrepid_example_Drivers_Example_09.exe deg NX NY verbose\n\n"; 00108 std::cout <<" where \n"; 00109 std::cout <<" int deg - polynomial degree to be used (assumed > 1) \n"; 00110 std::cout <<" int NX - num intervals in x direction (assumed box domain, 0,1) \n"; 00111 std::cout <<" int NY - 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 Quadrilateral 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 00152 00153 // *********************************** CELL TOPOLOGY ********************************** 00154 00155 // Get cell topology for base hexahedron 00156 typedef shards::CellTopology CellTopology; 00157 CellTopology quad_4(shards::getCellTopologyData<shards::Quadrilateral<4> >() ); 00158 00159 // Get dimensions 00160 int numNodesPerElem = quad_4.getNodeCount(); 00161 int spaceDim = quad_4.getDimension(); 00162 00163 // *********************************** GENERATE MESH ************************************ 00164 00165 *outStream << "Generating mesh ... \n\n"; 00166 00167 *outStream << " NX" << " NY\n"; 00168 *outStream << std::setw(5) << NX << 00169 std::setw(5) << NY << "\n\n"; 00170 00171 // Print mesh information 00172 int numElems = NX*NY; 00173 int numNodes = (NX+1)*(NY+1); 00174 *outStream << " Number of Elements: " << numElems << " \n"; 00175 *outStream << " Number of Nodes: " << numNodes << " \n\n"; 00176 00177 // Square 00178 double leftX = 0.0, rightX = 1.0; 00179 double leftY = 0.0, rightY = 1.0; 00180 00181 // Mesh spacing 00182 double hx = (rightX-leftX)/((double)NX); 00183 double hy = (rightY-leftY)/((double)NY); 00184 00185 // Get nodal coordinates 00186 FieldContainer<double> nodeCoord(numNodes, spaceDim); 00187 FieldContainer<int> nodeOnBoundary(numNodes); 00188 int inode = 0; 00189 for (int j=0; j<NY+1; j++) { 00190 for (int i=0; i<NX+1; i++) { 00191 nodeCoord(inode,0) = leftX + (double)i*hx; 00192 nodeCoord(inode,1) = leftY + (double)j*hy; 00193 if (j==0 || i==0 || j==NY || i==NX){ 00194 nodeOnBoundary(inode)=1; 00195 } 00196 else { 00197 nodeOnBoundary(inode)=0; 00198 } 00199 inode++; 00200 } 00201 } 00202 #define DUMP_DATA 00203 #ifdef DUMP_DATA 00204 // Print nodal coords 00205 ofstream fcoordout("coords.dat"); 00206 for (int i=0; i<numNodes; i++) { 00207 fcoordout << nodeCoord(i,0) <<" "; 00208 fcoordout << nodeCoord(i,1) <<"\n"; 00209 } 00210 fcoordout.close(); 00211 #endif 00212 00213 00214 // Element to Node map 00215 // We'll keep it around, but this is only the DOFMap if you are in the lowest order case. 00216 FieldContainer<int> elemToNode(numElems, numNodesPerElem); 00217 int ielem = 0; 00218 for (int j=0; j<NY; j++) { 00219 for (int i=0; i<NX; i++) { 00220 elemToNode(ielem,0) = (NX + 1)*j + i; 00221 elemToNode(ielem,1) = (NX + 1)*j + i + 1; 00222 elemToNode(ielem,2) = (NX + 1)*(j + 1) + i + 1; 00223 elemToNode(ielem,3) = (NX + 1)*(j + 1) + i; 00224 ielem++; 00225 } 00226 } 00227 #ifdef DUMP_DATA 00228 // Output connectivity 00229 ofstream fe2nout("elem2node.dat"); 00230 for (int j=0; j<NY; j++) { 00231 for (int i=0; i<NX; i++) { 00232 int ielem = i + j * NX; 00233 for (int m=0; m<numNodesPerElem; m++){ 00234 fe2nout << elemToNode(ielem,m) <<" "; 00235 } 00236 fe2nout <<"\n"; 00237 } 00238 } 00239 fe2nout.close(); 00240 #endif 00241 00242 00243 // ************************************ CUBATURE ************************************** 00244 *outStream << "Getting cubature ... \n\n"; 00245 00246 // Get numerical integration points and weights 00247 // I only need this on the line since I'm doing tensor products 00248 DefaultCubatureFactory<double> cubFactory; 00249 00250 Teuchos::RCP<Cubature<double,FieldContainer<double>,FieldContainer<double> > > glcub 00251 = Teuchos::rcp(new CubaturePolylib<double,FieldContainer<double>,FieldContainer<double> >(2*deg-1,PL_GAUSS_LOBATTO) ); 00252 00253 const int numCubPoints = glcub->getNumPoints(); 00254 00255 FieldContainer<double> cubPoints1D(numCubPoints, 1); 00256 FieldContainer<double> cubWeights1D(numCubPoints); 00257 00258 glcub->getCubature(cubPoints1D,cubWeights1D); 00259 00260 00261 // ************************************** BASIS *************************************** 00262 *outStream << "Getting basis ... \n\n"; 00263 00264 // Define basis: I only need this on the line also 00265 Basis_HGRAD_LINE_Cn_FEM<double, FieldContainer<double> > lineHGradBasis(deg,POINTTYPE_SPECTRAL); 00266 int numLineFieldsG = lineHGradBasis.getCardinality(); 00267 FieldContainer<double> lineGrads(numLineFieldsG, numCubPoints, 1); 00268 00269 // Evaluate basis values and gradients at cubature points 00270 lineHGradBasis.getValues(lineGrads, cubPoints1D, OPERATOR_GRAD); 00271 00272 // ************************************** LTG mapping ********************************** 00273 FieldContainer<int> ltgMapping(numElems,numLineFieldsG*numLineFieldsG); 00274 const int numDOF = (NX*deg+1)*(NY*deg+1); 00275 ielem=0; 00276 for (int j=0;j<NY;j++) { 00277 for (int i=0;i<NX;i++) { 00278 const int start = deg * j * ( NX * deg + 1 ) + i * deg; 00279 // loop over local dof on this cell 00280 int local_dof_cur=0; 00281 for (int vertical=0;vertical<=deg;vertical++) { 00282 for (int horizontal=0;horizontal<=deg;horizontal++) { 00283 ltgMapping(ielem,local_dof_cur) = start + vertical*(NX*deg+1)+horizontal; 00284 local_dof_cur++; 00285 } 00286 } 00287 ielem++; 00288 } 00289 } 00290 #ifdef DUMP_DATA 00291 // Output ltg mapping 00292 ofstream ltgout("ltg.dat"); 00293 for (int j=0; j<NY; j++) { 00294 for (int i=0; i<NX; i++) { 00295 int ielem = i + j * NX; 00296 for (int m=0; m<numLineFieldsG; m++){ 00297 ltgout << ltgMapping(ielem,m) <<" "; 00298 } 00299 ltgout <<"\n"; 00300 } 00301 } 00302 ltgout.close(); 00303 #endif 00304 00305 00306 // Global arrays in Epetra format 00307 Epetra_SerialComm Comm; 00308 Epetra_Map globalMapG(numDOF, 0, Comm); 00309 00310 Epetra_FEVector u(globalMapG); 00311 Epetra_FEVector Ku(globalMapG); 00312 00313 u.Random(); 00314 00315 00316 // ************************** Compute element HGrad stiffness matrices ******************************* 00317 // // Get vertices of all the cells 00318 // for (int i=0;i<numElems;i++) 00319 // { 00320 // for (int j=0;j<4;j++) 00321 // { 00322 // const int nodeCur = elemToNode(i,j); 00323 // for (int k=0;k<spaceDim;k++) 00324 // { 00325 // cellVertices(i,j,k) = nodeCoord(nodeCur,k); 00326 // } 00327 // } 00328 // } 00329 00330 FieldContainer<double> uScattered(numElems,numLineFieldsG*numLineFieldsG); 00331 FieldContainer<double> KuScattered(numElems,numLineFieldsG*numLineFieldsG); 00332 00333 // need storage for derivatives of u on each cell 00334 // the number of line dof should be the same as the 00335 // number of cub points. 00336 // This is indexed by Du(q2,q1) 00337 FieldContainer<double> Du(numCubPoints,numCubPoints); 00338 00339 00340 00341 double *uVals = u[0]; 00342 double *KuVals = Ku[0]; 00343 Epetra_Time scatterTime(Comm); 00344 *outStream << "Scattering\n"; 00345 // Scatter 00346 for (int k=0; k<numElems; k++) 00347 { 00348 for (int i=0;i<numLineFieldsG*numLineFieldsG;i++) 00349 { 00350 uScattered(k,i) = uVals[ltgMapping(k,i)]; 00351 } 00352 } 00353 const double scatTime = scatterTime.ElapsedTime(); 00354 *outStream << "Scattered in time " << scatTime << "\n"; 00355 00356 Epetra_Time applyTimer(Comm); 00357 00358 uScattered.resize(numElems,numLineFieldsG,numLineFieldsG); 00359 00360 for (int k=0;k<numElems;k++) 00361 { 00362 // local operation: x-derivative term of stiffness matrix 00363 // evaluate x derivative of u on cell k. 00364 for (int j=0;j<numLineFieldsG;j++) 00365 { 00366 for (int i=0;i<numLineFieldsG;i++) 00367 { 00368 Du(j,i) = 0.0; 00369 for (int q=0;q<numLineFieldsG;q++) 00370 { 00371 Du(j,i) += uScattered(k,j,i) * lineGrads(i,q,0); 00372 } 00373 } 00374 } 00375 00376 // initialize Ku 00377 for (int i=0;i<numLineFieldsG*numLineFieldsG;i++) 00378 { 00379 KuScattered(k,i) = 0.0; 00380 } 00381 00382 // loop over basis functions for x term 00383 int cur = 0; 00384 for (int j=0;j<numLineFieldsG;j++) 00385 { 00386 for (int i=0;i<numLineFieldsG;i++) 00387 { 00388 // do the quadrature 00389 for (int q1=0;q1<numCubPoints;q1++) 00390 { 00391 KuScattered(k,cur) += cubWeights1D(j) * cubWeights1D(q1) * Du(j,q1) * lineGrads(i,q1,0); 00392 } 00393 cur ++; 00394 } 00395 } 00396 00397 // local operation: y-derivative term of stiffness matrix, store in Du 00398 for (int j=0;j<numLineFieldsG;j++) 00399 { 00400 for (int i=0;i<numLineFieldsG;i++) 00401 { 00402 Du(j,i) = 0.0; 00403 for (int q=0;q<numLineFieldsG;q++) 00404 { 00405 Du(j,i) += uScattered(k,j,i) * lineGrads(j,q,0); 00406 } 00407 } 00408 } 00409 00410 00411 // evaluate y-derivatives of u 00412 cur = 0; 00413 for (int j=0;j<numLineFieldsG;j++) 00414 { 00415 for (int i=0;i<numLineFieldsG;i++) 00416 { 00417 for (int q2=0;q2<numCubPoints;q2++) 00418 { 00419 KuScattered(k,cur) += cubWeights1D(q2) * Du(q2,i) * lineGrads(j,q2,0) * cubWeights1D(i); 00420 } 00421 } 00422 } 00423 } 00424 00425 uScattered.resize(numElems,numLineFieldsG*numLineFieldsG); 00426 00427 const double applyTime = applyTimer.ElapsedTime(); 00428 00429 *outStream << "Local application: " << applyTime << "\n"; 00430 00431 // gather 00432 Epetra_Time gatherTimer(Comm); 00433 // Gather 00434 for (int k=0;k<numElems;k++) 00435 { 00436 for (int i=0;i<numLineFieldsG*numLineFieldsG;i++) 00437 { 00438 KuVals[ltgMapping(k,i)] += KuScattered(k,i); 00439 } 00440 } 00441 00442 const double gatherTime = gatherTimer.ElapsedTime(); 00443 *outStream << "Gathered in " << gatherTime << "\n"; 00444 00445 00446 00447 #ifdef DUMP_DATA 00448 // Dump matrices to disk 00449 // EpetraExt::RowMatrixToMatlabFile("stiff_matrix.dat",StiffMatrix); 00450 // EpetraExt::MultiVectorToMatrixMarketFile("rhs_vector.dat",rhs,0,0,false); 00451 #endif 00452 00453 00454 std::cout << "End Result: TEST PASSED\n"; 00455 00456 // reset format state of std::cout 00457 std::cout.copyfmt(oldFormatState); 00458 00459 return 0; 00460 } 00461
1.7.4