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