|
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.hpp" 00038 #include "Intrepid_DefaultCubatureFactory.hpp" 00039 #include "Intrepid_RealSpaceTools.hpp" 00040 #include "Intrepid_ArrayTools.hpp" 00041 #include "Intrepid_FunctionSpaceTools.hpp" 00042 #include "Intrepid_CellTools.hpp" 00043 #include "Teuchos_oblackholestream.hpp" 00044 #include "Teuchos_RCP.hpp" 00045 #include "Teuchos_GlobalMPISession.hpp" 00046 #include "Teuchos_SerialDenseMatrix.hpp" 00047 #include "Teuchos_SerialDenseVector.hpp" 00048 #include "Teuchos_LAPACK.hpp" 00049 00050 using namespace std; 00051 using namespace Intrepid; 00052 00053 void rhsFunc(FieldContainer<double> &, const FieldContainer<double> &, int, int, int); 00054 void neumann(FieldContainer<double> & , 00055 const FieldContainer<double> & , 00056 const FieldContainer<double> & , 00057 const shards::CellTopology & , 00058 int, int, int, int); 00059 void u_exact(FieldContainer<double> &, const FieldContainer<double> &, int, int, int); 00060 00062 void rhsFunc(FieldContainer<double> & result, 00063 const FieldContainer<double> & points, 00064 int xd, 00065 int yd, 00066 int zd) { 00067 00068 int x = 0, y = 1, z = 2; 00069 00070 // second x-derivatives of u 00071 if (xd > 1) { 00072 for (int cell=0; cell<result.dimension(0); cell++) { 00073 for (int pt=0; pt<result.dimension(1); pt++) { 00074 result(cell,pt) = - xd*(xd-1)*std::pow(points(cell,pt,x), xd-2) * 00075 std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd); 00076 } 00077 } 00078 } 00079 00080 // second y-derivatives of u 00081 if (yd > 1) { 00082 for (int cell=0; cell<result.dimension(0); cell++) { 00083 for (int pt=0; pt<result.dimension(1); pt++) { 00084 result(cell,pt) -= yd*(yd-1)*std::pow(points(cell,pt,y), yd-2) * 00085 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd); 00086 } 00087 } 00088 } 00089 00090 // second z-derivatives of u 00091 if (zd > 1) { 00092 for (int cell=0; cell<result.dimension(0); cell++) { 00093 for (int pt=0; pt<result.dimension(1); pt++) { 00094 result(cell,pt) -= zd*(zd-1)*std::pow(points(cell,pt,z), zd-2) * 00095 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd); 00096 } 00097 } 00098 } 00099 00100 // add u 00101 for (int cell=0; cell<result.dimension(0); cell++) { 00102 for (int pt=0; pt<result.dimension(1); pt++) { 00103 result(cell,pt) += std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd); 00104 } 00105 } 00106 00107 } 00108 00109 00111 void neumann(FieldContainer<double> & result, 00112 const FieldContainer<double> & points, 00113 const FieldContainer<double> & jacs, 00114 const shards::CellTopology & parentCell, 00115 int sideOrdinal, int xd, int yd, int zd) { 00116 00117 int x = 0, y = 1, z = 2; 00118 00119 int numCells = result.dimension(0); 00120 int numPoints = result.dimension(1); 00121 00122 FieldContainer<double> grad_u(numCells, numPoints, 3); 00123 FieldContainer<double> side_normals(numCells, numPoints, 3); 00124 FieldContainer<double> normal_lengths(numCells, numPoints); 00125 00126 // first x-derivatives of u 00127 if (xd > 0) { 00128 for (int cell=0; cell<numCells; cell++) { 00129 for (int pt=0; pt<numPoints; pt++) { 00130 grad_u(cell,pt,x) = xd*std::pow(points(cell,pt,x), xd-1) * 00131 std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd); 00132 } 00133 } 00134 } 00135 00136 // first y-derivatives of u 00137 if (yd > 0) { 00138 for (int cell=0; cell<numCells; cell++) { 00139 for (int pt=0; pt<numPoints; pt++) { 00140 grad_u(cell,pt,y) = yd*std::pow(points(cell,pt,y), yd-1) * 00141 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd); 00142 } 00143 } 00144 } 00145 00146 // first z-derivatives of u 00147 if (zd > 0) { 00148 for (int cell=0; cell<numCells; cell++) { 00149 for (int pt=0; pt<numPoints; pt++) { 00150 grad_u(cell,pt,z) = zd*std::pow(points(cell,pt,z), zd-1) * 00151 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd); 00152 } 00153 } 00154 } 00155 00156 CellTools<double>::getPhysicalSideNormals(side_normals, jacs, sideOrdinal, parentCell); 00157 00158 // scale normals 00159 RealSpaceTools<double>::vectorNorm(normal_lengths, side_normals, NORM_TWO); 00160 FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals, true); 00161 00162 FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals); 00163 00164 } 00165 00167 void u_exact(FieldContainer<double> & result, const FieldContainer<double> & points, int xd, int yd, int zd) { 00168 int x = 0, y = 1, z = 2; 00169 for (int cell=0; cell<result.dimension(0); cell++) { 00170 for (int pt=0; pt<result.dimension(1); pt++) { 00171 result(cell,pt) = std::pow(points(pt,x), xd)*std::pow(points(pt,y), yd)*std::pow(points(pt,z), zd); 00172 } 00173 } 00174 } 00175 00176 00177 00178 00179 int main(int argc, char *argv[]) { 00180 00181 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00182 00183 // This little trick lets us print to std::cout only if 00184 // a (dummy) command-line argument is provided. 00185 int iprint = argc - 1; 00186 Teuchos::RCP<std::ostream> outStream; 00187 Teuchos::oblackholestream bhs; // outputs nothing 00188 if (iprint > 0) 00189 outStream = Teuchos::rcp(&std::cout, false); 00190 else 00191 outStream = Teuchos::rcp(&bhs, false); 00192 00193 // Save the format state of the original std::cout. 00194 Teuchos::oblackholestream oldFormatState; 00195 oldFormatState.copyfmt(std::cout); 00196 00197 *outStream \ 00198 << "===============================================================================\n" \ 00199 << "| |\n" \ 00200 << "| Unit Test (Basis_HGRAD_TET_Cn_FEM) |\n" \ 00201 << "| |\n" \ 00202 << "| 1) Patch test involving mass and stiffness matrices, |\n" \ 00203 << "| for the Neumann problem on a tetrahedral patch |\n" \ 00204 << "| Omega with boundary Gamma. |\n" \ 00205 << "| |\n" \ 00206 << "| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \ 00207 << "| |\n" \ 00208 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00209 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00210 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00211 << "| |\n" \ 00212 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00213 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00214 << "| |\n" \ 00215 << "===============================================================================\n"\ 00216 << "| TEST 1: Patch test |\n"\ 00217 << "===============================================================================\n"; 00218 00219 00220 int errorFlag = 0; 00221 00222 outStream -> precision(16); 00223 00224 00225 try { 00226 00227 int max_order = 5; // max total order of polynomial solution 00228 DefaultCubatureFactory<double> cubFactory; // create factory 00229 shards::CellTopology cell(shards::getCellTopologyData< shards::Tetrahedron<> >()); // create parent cell topology 00230 shards::CellTopology side(shards::getCellTopologyData< shards::Triangle<> >()); // create relevant subcell (side) topology 00231 int cellDim = cell.getDimension(); 00232 int sideDim = side.getDimension(); 00233 00234 // Define array containing points at which the solution is evaluated, on the reference tet. 00235 int numIntervals = 10; 00236 int numInterpPoints = ((numIntervals + 1)*(numIntervals + 2)*(numIntervals + 3))/6; 00237 FieldContainer<double> interp_points_ref(numInterpPoints, 3); 00238 int counter = 0; 00239 for (int k=0; k<=numIntervals; k++) { 00240 for (int j=0; j<=numIntervals; j++) { 00241 for (int i=0; i<=numIntervals; i++) { 00242 if (i+j+k <= numIntervals) { 00243 interp_points_ref(counter,0) = i*(1.0/numIntervals); 00244 interp_points_ref(counter,1) = j*(1.0/numIntervals); 00245 interp_points_ref(counter,2) = k*(1.0/numIntervals); 00246 counter++; 00247 } 00248 } 00249 } 00250 } 00251 00252 /* Definition of parent cell. */ 00253 FieldContainer<double> cell_nodes(1, 4, cellDim); 00254 // funky tet 00255 cell_nodes(0, 0, 0) = -1.0; 00256 cell_nodes(0, 0, 1) = -2.0; 00257 cell_nodes(0, 0, 2) = 0.0; 00258 cell_nodes(0, 1, 0) = 6.0; 00259 cell_nodes(0, 1, 1) = 2.0; 00260 cell_nodes(0, 1, 2) = 0.0; 00261 cell_nodes(0, 2, 0) = -5.0; 00262 cell_nodes(0, 2, 1) = 1.0; 00263 cell_nodes(0, 2, 2) = 0.0; 00264 cell_nodes(0, 3, 0) = -4.0; 00265 cell_nodes(0, 3, 1) = -1.0; 00266 cell_nodes(0, 3, 2) = 3.0; 00267 // perturbed reference tet 00268 /*cell_nodes(0, 0, 0) = 0.1; 00269 cell_nodes(0, 0, 1) = -0.1; 00270 cell_nodes(0, 0, 2) = 0.2; 00271 cell_nodes(0, 1, 0) = 1.2; 00272 cell_nodes(0, 1, 1) = -0.1; 00273 cell_nodes(0, 1, 2) = 0.05; 00274 cell_nodes(0, 2, 0) = 0.0; 00275 cell_nodes(0, 2, 1) = 0.9; 00276 cell_nodes(0, 2, 2) = 0.1; 00277 cell_nodes(0, 3, 0) = 0.1; 00278 cell_nodes(0, 3, 1) = -0.1; 00279 cell_nodes(0, 3, 2) = 1.1;*/ 00280 // reference tet 00281 /*cell_nodes(0, 0, 0) = 0.0; 00282 cell_nodes(0, 0, 1) = 0.0; 00283 cell_nodes(0, 0, 2) = 0.0; 00284 cell_nodes(0, 1, 0) = 1.0; 00285 cell_nodes(0, 1, 1) = 0.0; 00286 cell_nodes(0, 1, 2) = 0.0; 00287 cell_nodes(0, 2, 0) = 0.0; 00288 cell_nodes(0, 2, 1) = 1.0; 00289 cell_nodes(0, 2, 2) = 0.0; 00290 cell_nodes(0, 3, 0) = 0.0; 00291 cell_nodes(0, 3, 1) = 0.0; 00292 cell_nodes(0, 3, 2) = 1.0;*/ 00293 00294 FieldContainer<double> interp_points(1, numInterpPoints, cellDim); 00295 CellTools<double>::mapToPhysicalFrame(interp_points, interp_points_ref, cell_nodes, cell); 00296 interp_points.resize(numInterpPoints, cellDim); 00297 00298 // we test two types of bases 00299 EPointType pointtype[] = {POINTTYPE_EQUISPACED, POINTTYPE_WARPBLEND}; 00300 for (int ptype=0; ptype < 2; ptype++) { 00301 00302 *outStream << "\nTesting bases with " << EPointTypeToString(pointtype[ptype]) << ":\n"; 00303 00304 for (int x_order=0; x_order <= max_order; x_order++) { 00305 for (int y_order=0; y_order <= max_order-x_order; y_order++) { 00306 for (int z_order=0; z_order <= max_order-x_order-y_order; z_order++) { 00307 00308 // evaluate exact solution 00309 FieldContainer<double> exact_solution(1, numInterpPoints); 00310 u_exact(exact_solution, interp_points, x_order, y_order, z_order); 00311 00312 int total_order = std::max(x_order + y_order + z_order, 1); 00313 00314 for (int basis_order=total_order; basis_order <= max_order; basis_order++) { 00315 00316 // set test tolerance; 00317 double zero = basis_order*basis_order*basis_order*100*INTREPID_TOL; 00318 00319 //create basis 00320 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis = 00321 Teuchos::rcp(new Basis_HGRAD_TET_Cn_FEM<double,FieldContainer<double> >(basis_order, pointtype[ptype]) ); 00322 int numFields = basis->getCardinality(); 00323 00324 // create cubatures 00325 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*basis_order); 00326 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*basis_order); 00327 int numCubPointsCell = cellCub->getNumPoints(); 00328 int numCubPointsSide = sideCub->getNumPoints(); 00329 00330 /* Computational arrays. */ 00331 /* Section 1: Related to parent cell integration. */ 00332 FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim); 00333 FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim); 00334 FieldContainer<double> cub_weights_cell(numCubPointsCell); 00335 FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim); 00336 FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim); 00337 FieldContainer<double> jacobian_det_cell(1, numCubPointsCell); 00338 FieldContainer<double> weighted_measure_cell(1, numCubPointsCell); 00339 00340 FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell); 00341 FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell); 00342 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell); 00343 FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim); 00344 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim); 00345 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim); 00346 FieldContainer<double> fe_matrix(1, numFields, numFields); 00347 00348 FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell); 00349 FieldContainer<double> rhs_and_soln_vector(1, numFields); 00350 00351 /* Section 2: Related to subcell (side) integration. */ 00352 unsigned numSides = 4; 00353 FieldContainer<double> cub_points_side(numCubPointsSide, sideDim); 00354 FieldContainer<double> cub_weights_side(numCubPointsSide); 00355 FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim); 00356 FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim); 00357 FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim); 00358 FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide); 00359 FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide); 00360 00361 FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide); 00362 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide); 00363 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide); 00364 FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide); 00365 FieldContainer<double> neumann_fields_per_side(1, numFields); 00366 00367 /* Section 3: Related to global interpolant. */ 00368 FieldContainer<double> value_of_basis_at_interp_points_ref(numFields, numInterpPoints); 00369 FieldContainer<double> transformed_value_of_basis_at_interp_points_ref(1, numFields, numInterpPoints); 00370 FieldContainer<double> interpolant(1, numInterpPoints); 00371 00372 FieldContainer<int> ipiv(numFields); 00373 00374 00375 00376 /******************* START COMPUTATION ***********************/ 00377 00378 // get cubature points and weights 00379 cellCub->getCubature(cub_points_cell, cub_weights_cell); 00380 00381 // compute geometric cell information 00382 CellTools<double>::setJacobian(jacobian_cell, cub_points_cell, cell_nodes, cell); 00383 CellTools<double>::setJacobianInv(jacobian_inv_cell, jacobian_cell); 00384 CellTools<double>::setJacobianDet(jacobian_det_cell, jacobian_cell); 00385 00386 // compute weighted measure 00387 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell); 00388 00390 // Computing mass matrices: 00391 // tabulate values of basis functions at (reference) cubature points 00392 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE); 00393 00394 // transform values of basis functions 00395 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell, 00396 value_of_basis_at_cub_points_cell); 00397 00398 // multiply with weighted measure 00399 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell, 00400 weighted_measure_cell, 00401 transformed_value_of_basis_at_cub_points_cell); 00402 00403 // compute mass matrices 00404 FunctionSpaceTools::integrate<double>(fe_matrix, 00405 transformed_value_of_basis_at_cub_points_cell, 00406 weighted_transformed_value_of_basis_at_cub_points_cell, 00407 COMP_BLAS); 00409 00411 // Computing stiffness matrices: 00412 // tabulate gradients of basis functions at (reference) cubature points 00413 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD); 00414 00415 // transform gradients of basis functions 00416 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell, 00417 jacobian_inv_cell, 00418 grad_of_basis_at_cub_points_cell); 00419 00420 // multiply with weighted measure 00421 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell, 00422 weighted_measure_cell, 00423 transformed_grad_of_basis_at_cub_points_cell); 00424 00425 // compute stiffness matrices and sum into fe_matrix 00426 FunctionSpaceTools::integrate<double>(fe_matrix, 00427 transformed_grad_of_basis_at_cub_points_cell, 00428 weighted_transformed_grad_of_basis_at_cub_points_cell, 00429 COMP_BLAS, 00430 true); 00432 00434 // Computing RHS contributions: 00435 // map cell (reference) cubature points to physical space 00436 CellTools<double>::mapToPhysicalFrame(cub_points_cell_physical, cub_points_cell, cell_nodes, cell); 00437 00438 // evaluate rhs function 00439 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order, z_order); 00440 00441 // compute rhs 00442 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector, 00443 rhs_at_cub_points_cell_physical, 00444 weighted_transformed_value_of_basis_at_cub_points_cell, 00445 COMP_BLAS); 00446 00447 // compute neumann b.c. contributions and adjust rhs 00448 sideCub->getCubature(cub_points_side, cub_weights_side); 00449 for (unsigned i=0; i<numSides; i++) { 00450 // compute geometric cell information 00451 CellTools<double>::mapToReferenceSubcell(cub_points_side_refcell, cub_points_side, sideDim, (int)i, cell); 00452 CellTools<double>::setJacobian(jacobian_side_refcell, cub_points_side_refcell, cell_nodes, cell); 00453 CellTools<double>::setJacobianDet(jacobian_det_side_refcell, jacobian_side_refcell); 00454 00455 // compute weighted face measure 00456 FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_side_refcell, 00457 jacobian_side_refcell, 00458 cub_weights_side, 00459 i, 00460 cell); 00461 00462 // tabulate values of basis functions at side cubature points, in the reference parent cell domain 00463 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE); 00464 // transform 00465 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell, 00466 value_of_basis_at_cub_points_side_refcell); 00467 00468 // multiply with weighted measure 00469 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell, 00470 weighted_measure_side_refcell, 00471 transformed_value_of_basis_at_cub_points_side_refcell); 00472 00473 // compute Neumann data 00474 // map side cubature points in reference parent cell domain to physical space 00475 CellTools<double>::mapToPhysicalFrame(cub_points_side_physical, cub_points_side_refcell, cell_nodes, cell); 00476 // now compute data 00477 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell, 00478 cell, (int)i, x_order, y_order, z_order); 00479 00480 FunctionSpaceTools::integrate<double>(neumann_fields_per_side, 00481 neumann_data_at_cub_points_side_physical, 00482 weighted_transformed_value_of_basis_at_cub_points_side_refcell, 00483 COMP_BLAS); 00484 00485 // adjust RHS 00486 RealSpaceTools<double>::add(rhs_and_soln_vector, neumann_fields_per_side);; 00487 } 00489 00491 // Solution of linear system: 00492 int info = 0; 00493 Teuchos::LAPACK<int, double> solver; 00494 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info); 00496 00498 // Building interpolant: 00499 // evaluate basis at interpolation points 00500 basis->getValues(value_of_basis_at_interp_points_ref, interp_points_ref, OPERATOR_VALUE); 00501 // transform values of basis functions 00502 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points_ref, 00503 value_of_basis_at_interp_points_ref); 00504 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points_ref); 00506 00507 /******************* END COMPUTATION ***********************/ 00508 00509 RealSpaceTools<double>::subtract(interpolant, exact_solution); 00510 00511 *outStream << "\nRelative norm-2 error between exact solution polynomial of order (" 00512 << x_order << ", " << y_order << ", " << z_order 00513 << ") and finite element interpolant of order " << basis_order << ": " 00514 << RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) / 00515 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) << "\n"; 00516 00517 if (RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) / 00518 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) > zero) { 00519 *outStream << "\n\nPatch test failed for solution polynomial order (" 00520 << x_order << ", " << y_order << ", " << z_order << ") and basis order " << basis_order << "\n\n"; 00521 errorFlag++; 00522 } 00523 } // end for basis_order 00524 } // end for z_order 00525 } // end for y_order 00526 } // end for x_order 00527 } // end for ptype 00528 00529 } 00530 // Catch unexpected errors 00531 catch (std::logic_error err) { 00532 *outStream << err.what() << "\n\n"; 00533 errorFlag = -1000; 00534 }; 00535 00536 if (errorFlag != 0) 00537 std::cout << "End Result: TEST FAILED\n"; 00538 else 00539 std::cout << "End Result: TEST PASSED\n"; 00540 00541 // reset format state of std::cout 00542 std::cout.copyfmt(oldFormatState); 00543 00544 return errorFlag; 00545 }
1.7.4