|
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 00036 #include "Intrepid_FieldContainer.hpp" 00037 #include "Intrepid_HGRAD_HEX_Cn_FEM.hpp" 00038 #include "Intrepid_HDIV_HEX_In_FEM.hpp" 00039 #include "Intrepid_DefaultCubatureFactory.hpp" 00040 #include "Intrepid_PointTools.hpp" 00041 #include "Intrepid_RealSpaceTools.hpp" 00042 #include "Intrepid_ArrayTools.hpp" 00043 #include "Intrepid_FunctionSpaceTools.hpp" 00044 #include "Intrepid_CellTools.hpp" 00045 #include "Teuchos_oblackholestream.hpp" 00046 #include "Teuchos_RCP.hpp" 00047 #include "Teuchos_GlobalMPISession.hpp" 00048 #include "Teuchos_SerialDenseMatrix.hpp" 00049 #include "Teuchos_SerialDenseVector.hpp" 00050 #include "Teuchos_LAPACK.hpp" 00051 00052 using namespace std; 00053 using namespace Intrepid; 00054 00055 void rhsFunc( FieldContainer<double> &, const FieldContainer<double> &, int, int, int ); 00056 void u_exact( FieldContainer<double> &, const FieldContainer<double> &, int, int, int ); 00057 00058 // This is the rhs for (div tau,w) = (f,w), 00059 // which makes f the negative Laplacian of scalar solution 00060 void rhsFunc( FieldContainer<double> &result, 00061 const FieldContainer<double> &points, 00062 int xd, 00063 int yd , 00064 int zd) 00065 { 00066 for (int cell=0;cell<result.dimension(0);cell++) { 00067 for (int pt=0;pt<result.dimension(1);pt++) { 00068 result(cell,pt) = 0.0; 00069 if (xd >=2) { 00070 result(cell,pt) += xd*(xd-1)*pow(points(cell,pt,0),xd-2)*pow(points(cell,pt,1),yd) 00071 *pow(points(cell,pt,2),zd); 00072 } 00073 if (yd >=2) { 00074 result(cell,pt) += yd*(yd-1)*pow(points(cell,pt,0),xd)*pow(points(cell,pt,1),yd-2) 00075 *pow(points(cell,pt,2),zd); 00076 } 00077 if (zd>=2) { 00078 result(cell,pt) += zd*(zd-1)*pow(points(cell,pt,0),xd)*pow(points(cell,pt,1),yd) 00079 *pow(points(cell,pt,2),zd-2); 00080 00081 } 00082 } 00083 } 00084 } 00085 00086 void u_exact( FieldContainer<double> &result, 00087 const FieldContainer<double> &points, 00088 int xd, 00089 int yd, 00090 int zd) 00091 { 00092 for (int cell=0;cell<result.dimension(0);cell++){ 00093 for (int pt=0;pt<result.dimension(1);pt++) { 00094 result(cell,pt) = std::pow(points(cell,pt,0),xd)*std::pow(points(cell,pt,1),yd) 00095 *std::pow(points(cell,pt,2),zd); 00096 } 00097 } 00098 return; 00099 } 00100 00101 int main(int argc, char *argv[]) { 00102 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00103 00104 // This little trick lets us print to std::cout only if 00105 // a (dummy) command-line argument is provided. 00106 int iprint = argc - 1; 00107 Teuchos::RCP<std::ostream> outStream; 00108 Teuchos::oblackholestream bhs; // outputs nothing 00109 if (iprint > 0) 00110 outStream = Teuchos::rcp(&std::cout, false); 00111 else 00112 outStream = Teuchos::rcp(&bhs, false); 00113 00114 // Save the format state of the original std::cout. 00115 Teuchos::oblackholestream oldFormatState; 00116 oldFormatState.copyfmt(std::cout); 00117 00118 *outStream \ 00119 << "===============================================================================\n" \ 00120 << "| |\n" \ 00121 << "| Unit Test (Basis_HGRAD_HEX_In_FEM) |\n" \ 00122 << "| |\n" \ 00123 << "| 1) Patch test involving H(div) matrices |\n" \ 00124 << "| for the Dirichlet problem on a hexahedron |\n" \ 00125 << "| Omega with boundary Gamma. |\n" \ 00126 << "| |\n" \ 00127 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00128 << "| Robert Kirby (robert.c.kirby@ttu.edu), |\n" \ 00129 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00130 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00131 << "| |\n" \ 00132 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00133 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00134 << "| |\n" \ 00135 << "===============================================================================\n" \ 00136 << "| TEST 1: Patch test |\n" \ 00137 << "===============================================================================\n"; 00138 00139 00140 int errorFlag = 0; 00141 00142 outStream -> precision(16); 00143 00144 try { 00145 DefaultCubatureFactory<double> cubFactory; // create cubature factory 00146 shards::CellTopology cell(shards::getCellTopologyData< shards::Hexahedron<> >()); // create parent cell topology 00147 shards::CellTopology side(shards::getCellTopologyData< shards::Quadrilateral<> >()); // create relevant subcell (side) topology 00148 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >() ); // for getting points to construct the basis 00149 00150 int cellDim = cell.getDimension(); 00151 int sideDim = side.getDimension(); 00152 00153 int min_order = 0; 00154 int max_order = 3; 00155 00156 int numIntervals = 2; 00157 int numInterpPoints = (numIntervals + 1)*(numIntervals + 1)*(numIntervals+1); 00158 FieldContainer<double> interp_points_ref(numInterpPoints, cellDim); 00159 int counter = 0; 00160 for (int k=0;k<numIntervals;k++) { 00161 for (int j=0; j<=numIntervals; j++) { 00162 for (int i=0; i<=numIntervals; i++) { 00163 interp_points_ref(counter,0) = i*(1.0/numIntervals); 00164 interp_points_ref(counter,1) = j*(1.0/numIntervals); 00165 interp_points_ref(counter,2) = k*(1.0/numIntervals); 00166 counter++; 00167 } 00168 } 00169 } 00170 00171 for (int basis_order=min_order;basis_order<=max_order;basis_order++) { 00172 // create bases 00173 // get the points for the vector basis 00174 Teuchos::RCP<Basis<double,FieldContainer<double> > > vectorBasis = 00175 Teuchos::rcp(new Basis_HDIV_HEX_In_FEM<double,FieldContainer<double> >(basis_order+1,POINTTYPE_SPECTRAL) ); 00176 00177 Teuchos::RCP<Basis<double,FieldContainer<double> > > scalarBasis = 00178 Teuchos::rcp(new Basis_HGRAD_HEX_Cn_FEM<double,FieldContainer<double> >(basis_order,POINTTYPE_SPECTRAL) ); 00179 00180 int numVectorFields = vectorBasis->getCardinality(); 00181 int numScalarFields = scalarBasis->getCardinality(); 00182 int numTotalFields = numVectorFields + numScalarFields; 00183 00184 // create cubatures 00185 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*(basis_order+1)); 00186 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*(basis_order+1)); 00187 00188 int numCubPointsCell = cellCub->getNumPoints(); 00189 int numCubPointsSide = sideCub->getNumPoints(); 00190 00191 // hold cubature information 00192 FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim); 00193 FieldContainer<double> cub_weights_cell(numCubPointsCell); 00194 FieldContainer<double> cub_points_side( numCubPointsSide, sideDim ); 00195 FieldContainer<double> cub_weights_side( numCubPointsSide ); 00196 FieldContainer<double> cub_points_side_refcell( numCubPointsSide , cellDim ); 00197 00198 // hold basis function information on refcell 00199 FieldContainer<double> value_of_v_basis_at_cub_points_cell(numVectorFields, numCubPointsCell, cellDim ); 00200 FieldContainer<double> w_value_of_v_basis_at_cub_points_cell(1, numVectorFields, numCubPointsCell, cellDim); 00201 FieldContainer<double> div_of_v_basis_at_cub_points_cell( numVectorFields, numCubPointsCell ); 00202 FieldContainer<double> w_div_of_v_basis_at_cub_points_cell( 1, numVectorFields , numCubPointsCell ); 00203 FieldContainer<double> value_of_s_basis_at_cub_points_cell(numScalarFields,numCubPointsCell); 00204 FieldContainer<double> w_value_of_s_basis_at_cub_points_cell(1,numScalarFields,numCubPointsCell); 00205 00206 // containers for side integration: 00207 // I just need the normal component of the vector basis 00208 // and the exact solution at the cub points 00209 FieldContainer<double> value_of_v_basis_at_cub_points_side(numVectorFields,numCubPointsSide,cellDim); 00210 FieldContainer<double> n_of_v_basis_at_cub_points_side(numVectorFields,numCubPointsSide); 00211 FieldContainer<double> w_n_of_v_basis_at_cub_points_side(1,numVectorFields,numCubPointsSide); 00212 FieldContainer<double> diri_data_at_cub_points_side(1,numCubPointsSide); 00213 FieldContainer<double> side_normal(cellDim); 00214 00215 // holds rhs data 00216 FieldContainer<double> rhs_at_cub_points_cell(1,numCubPointsCell); 00217 00218 // FEM matrices and vectors 00219 FieldContainer<double> fe_matrix_M(1,numVectorFields,numVectorFields); 00220 FieldContainer<double> fe_matrix_B(1,numVectorFields,numScalarFields); 00221 FieldContainer<double> fe_matrix(1,numTotalFields,numTotalFields); 00222 00223 FieldContainer<double> rhs_vector_vec(1,numVectorFields); 00224 FieldContainer<double> rhs_vector_scal(1,numScalarFields); 00225 FieldContainer<double> rhs_and_soln_vec(1,numTotalFields); 00226 00227 FieldContainer<int> ipiv(numTotalFields); 00228 FieldContainer<double> value_of_s_basis_at_interp_points( numScalarFields , numInterpPoints); 00229 FieldContainer<double> interpolant( 1 , numInterpPoints ); 00230 00231 // set test tolerance 00232 double zero = (basis_order+1)*(basis_order+1)*1000.0*INTREPID_TOL; 00233 00234 // build matrices outside the loop, and then just do the rhs 00235 // for each iteration 00236 00237 cellCub->getCubature(cub_points_cell, cub_weights_cell); 00238 sideCub->getCubature(cub_points_side, cub_weights_side); 00239 00240 // need the vector basis & its divergences 00241 vectorBasis->getValues(value_of_v_basis_at_cub_points_cell, 00242 cub_points_cell, 00243 OPERATOR_VALUE); 00244 vectorBasis->getValues(div_of_v_basis_at_cub_points_cell, 00245 cub_points_cell, 00246 OPERATOR_DIV); 00247 00248 // need the scalar basis as well 00249 scalarBasis->getValues(value_of_s_basis_at_cub_points_cell, 00250 cub_points_cell, 00251 OPERATOR_VALUE); 00252 00253 // construct mass matrix 00254 cub_weights_cell.resize(1,numCubPointsCell); 00255 FunctionSpaceTools::multiplyMeasure<double>(w_value_of_v_basis_at_cub_points_cell , 00256 cub_weights_cell , 00257 value_of_v_basis_at_cub_points_cell ); 00258 cub_weights_cell.resize(numCubPointsCell); 00259 00260 00261 value_of_v_basis_at_cub_points_cell.resize( 1 , numVectorFields , numCubPointsCell , cellDim ); 00262 FunctionSpaceTools::integrate<double>(fe_matrix_M, 00263 w_value_of_v_basis_at_cub_points_cell , 00264 value_of_v_basis_at_cub_points_cell , 00265 COMP_BLAS ); 00266 value_of_v_basis_at_cub_points_cell.resize( numVectorFields , numCubPointsCell , cellDim ); 00267 00268 // div matrix 00269 cub_weights_cell.resize(1,numCubPointsCell); 00270 FunctionSpaceTools::multiplyMeasure<double>(w_div_of_v_basis_at_cub_points_cell, 00271 cub_weights_cell, 00272 div_of_v_basis_at_cub_points_cell); 00273 cub_weights_cell.resize(numCubPointsCell); 00274 00275 value_of_s_basis_at_cub_points_cell.resize(1,numScalarFields,numCubPointsCell); 00276 FunctionSpaceTools::integrate<double>(fe_matrix_B, 00277 w_div_of_v_basis_at_cub_points_cell , 00278 value_of_s_basis_at_cub_points_cell , 00279 COMP_BLAS ); 00280 value_of_s_basis_at_cub_points_cell.resize(numScalarFields,numCubPointsCell); 00281 00282 for (int x_order=0;x_order<=basis_order;x_order++) { 00283 for (int y_order=0;y_order<=basis_order;y_order++) { 00284 for (int z_order=0;z_order<=basis_order;z_order++) { 00285 00286 00287 // reset global matrix since I destroyed it in LU factorization. 00288 fe_matrix.initialize(); 00289 // insert mass matrix into global matrix 00290 for (int i=0;i<numVectorFields;i++) { 00291 for (int j=0;j<numVectorFields;j++) { 00292 fe_matrix(0,i,j) = fe_matrix_M(0,i,j); 00293 } 00294 } 00295 00296 // insert div matrix into global matrix 00297 for (int i=0;i<numVectorFields;i++) { 00298 for (int j=0;j<numScalarFields;j++) { 00299 fe_matrix(0,i,numVectorFields+j)=-fe_matrix_B(0,i,j); 00300 fe_matrix(0,j+numVectorFields,i)=fe_matrix_B(0,i,j); 00301 } 00302 } 00303 00304 // clear old vector data 00305 rhs_vector_vec.initialize(); 00306 rhs_vector_scal.initialize(); 00307 rhs_and_soln_vec.initialize(); 00308 00309 // now get rhs vector 00310 // rhs_vector_scal is just (rhs,w) for w in the scalar basis 00311 // I already have the scalar basis tabulated. 00312 cub_points_cell.resize(1,numCubPointsCell,cellDim); 00313 rhsFunc(rhs_at_cub_points_cell, 00314 cub_points_cell, 00315 x_order, 00316 y_order, 00317 z_order); 00318 00319 cub_points_cell.resize(numCubPointsCell,cellDim); 00320 00321 cub_weights_cell.resize(1,numCubPointsCell); 00322 FunctionSpaceTools::multiplyMeasure<double>(w_value_of_s_basis_at_cub_points_cell, 00323 cub_weights_cell, 00324 value_of_s_basis_at_cub_points_cell); 00325 cub_weights_cell.resize(numCubPointsCell); 00326 FunctionSpaceTools::integrate<double>(rhs_vector_scal, 00327 rhs_at_cub_points_cell, 00328 w_value_of_s_basis_at_cub_points_cell, 00329 COMP_BLAS); 00330 00331 for (int i=0;i<numScalarFields;i++) { 00332 rhs_and_soln_vec(0,numVectorFields+i) = rhs_vector_scal(0,i); 00333 } 00334 00335 00336 // now get <u,v.n> on boundary 00337 for (unsigned side_cur=0;side_cur<6;side_cur++) { 00338 // map side cubature to current side 00339 CellTools<double>::mapToReferenceSubcell( cub_points_side_refcell , 00340 cub_points_side , 00341 sideDim , 00342 (int)side_cur , 00343 cell ); 00344 // Evaluate dirichlet data 00345 cub_points_side_refcell.resize(1,numCubPointsSide,cellDim); 00346 u_exact(diri_data_at_cub_points_side, 00347 cub_points_side_refcell,x_order,y_order,z_order); 00348 00349 cub_points_side_refcell.resize(numCubPointsSide,cellDim); 00350 00351 // get normal direction, this has the edge weight factored into it already 00352 CellTools<double>::getReferenceSideNormal(side_normal , 00353 (int)side_cur,cell ); 00354 00355 // v.n at cub points on side 00356 vectorBasis->getValues(value_of_v_basis_at_cub_points_side , 00357 cub_points_side_refcell , 00358 OPERATOR_VALUE ); 00359 00360 for (int i=0;i<numVectorFields;i++) { 00361 for (int j=0;j<numCubPointsSide;j++) { 00362 n_of_v_basis_at_cub_points_side(i,j) = 0.0; 00363 for (int k=0;k<cellDim;k++) { 00364 n_of_v_basis_at_cub_points_side(i,j) += side_normal(k) * 00365 value_of_v_basis_at_cub_points_side(i,j,k); 00366 } 00367 } 00368 } 00369 00370 cub_weights_side.resize(1,numCubPointsSide); 00371 FunctionSpaceTools::multiplyMeasure<double>(w_n_of_v_basis_at_cub_points_side, 00372 cub_weights_side, 00373 n_of_v_basis_at_cub_points_side); 00374 cub_weights_side.resize(numCubPointsSide); 00375 00376 FunctionSpaceTools::integrate<double>(rhs_vector_vec, 00377 diri_data_at_cub_points_side, 00378 w_n_of_v_basis_at_cub_points_side, 00379 COMP_BLAS, 00380 false); 00381 00382 for (int i=0;i<numVectorFields;i++) { 00383 rhs_and_soln_vec(0,i) -= rhs_vector_vec(0,i); 00384 } 00385 00386 } 00387 00388 // solve linear system 00389 int info = 0; 00390 Teuchos::LAPACK<int, double> solver; 00391 solver.GESV(numTotalFields, 1, &fe_matrix[0], numTotalFields, &ipiv(0), &rhs_and_soln_vec[0], 00392 numTotalFields, &info); 00393 00394 // compute interpolant; the scalar entries are last 00395 scalarBasis->getValues(value_of_s_basis_at_interp_points, 00396 interp_points_ref, 00397 OPERATOR_VALUE); 00398 for (int pt=0;pt<numInterpPoints;pt++) { 00399 interpolant(0,pt)=0.0; 00400 for (int i=0;i<numScalarFields;i++) { 00401 interpolant(0,pt) += rhs_and_soln_vec(0,numVectorFields+i) 00402 * value_of_s_basis_at_interp_points(i,pt); 00403 } 00404 } 00405 00406 interp_points_ref.resize(1,numInterpPoints,cellDim); 00407 // get exact solution for comparison 00408 FieldContainer<double> exact_solution(1,numInterpPoints); 00409 u_exact( exact_solution , interp_points_ref , x_order, y_order, z_order); 00410 interp_points_ref.resize(numInterpPoints,cellDim); 00411 00412 RealSpaceTools<double>::add(interpolant,exact_solution); 00413 00414 double nrm= RealSpaceTools<double>::vectorNorm(&interpolant[0],interpolant.dimension(1), NORM_TWO); 00415 00416 *outStream << "\nNorm-2 error between scalar components of exact solution of order (" 00417 << x_order << ", " << y_order << ", " << z_order 00418 << ") and finite element interpolant of order " << basis_order << ": " 00419 << nrm << "\n"; 00420 00421 if (nrm > zero) { 00422 *outStream << "\n\nPatch test failed for solution polynomial order (" 00423 << x_order << ", " << y_order << ", " << z_order << ") and basis order (scalar, vector) (" 00424 << basis_order << ", " << basis_order+1 << ")\n\n"; 00425 errorFlag++; 00426 } 00427 00428 } 00429 } 00430 } 00431 } 00432 00433 } 00434 00435 catch (std::logic_error err) { 00436 *outStream << err.what() << "\n\n"; 00437 errorFlag = -1000; 00438 }; 00439 00440 if (errorFlag != 0) 00441 std::cout << "End Result: TEST FAILED\n"; 00442 else 00443 std::cout << "End Result: TEST PASSED\n"; 00444 00445 // reset format state of std::cout 00446 std::cout.copyfmt(oldFormatState); 00447 00448 return errorFlag; 00449 }
1.7.4