00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00036 #include "Intrepid_FieldContainer.hpp"
00037 #include "Intrepid_HGRAD_TET_Cn_FEM_ORTH.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
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
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
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
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
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
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
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
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
00184
00185 int iprint = argc - 1;
00186 Teuchos::RCP<std::ostream> outStream;
00187 Teuchos::oblackholestream bhs;
00188 if (iprint > 0)
00189 outStream = Teuchos::rcp(&std::cout, false);
00190 else
00191 outStream = Teuchos::rcp(&bhs, false);
00192
00193
00194 Teuchos::oblackholestream oldFormatState;
00195 oldFormatState.copyfmt(std::cout);
00196
00197 *outStream \
00198 << "===============================================================================\n" \
00199 << "| |\n" \
00200 << "| Unit Test (Basis_HGRAD_TET_Cn_FEM_ORTH) |\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 = 7;
00228 DefaultCubatureFactory<double> cubFactory;
00229 shards::CellTopology cell(shards::getCellTopologyData< shards::Tetrahedron<> >());
00230 shards::CellTopology side(shards::getCellTopologyData< shards::Triangle<> >());
00231 int cellDim = cell.getDimension();
00232 int sideDim = side.getDimension();
00233
00234
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
00253 FieldContainer<double> cell_nodes(1, 4, cellDim);
00254
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
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
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 for (int x_order=0; x_order <= max_order; x_order++) {
00299 for (int y_order=0; y_order <= max_order-x_order; y_order++) {
00300 for (int z_order=0; z_order <= max_order-x_order-y_order; z_order++) {
00301
00302
00303 FieldContainer<double> exact_solution(1, numInterpPoints);
00304 u_exact(exact_solution, interp_points, x_order, y_order, z_order);
00305
00306 int total_order = std::max(x_order + y_order + z_order, 1);
00307
00308 for (int basis_order=total_order; basis_order <= max_order; basis_order++) {
00309
00310
00311 double zero = basis_order*basis_order*basis_order*100*INTREPID_TOL;
00312
00313
00314 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
00315 Teuchos::rcp(new Basis_HGRAD_TET_Cn_FEM_ORTH<double,FieldContainer<double> >(basis_order) );
00316 int numFields = basis->getCardinality();
00317
00318
00319 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*basis_order);
00320 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*basis_order);
00321 int numCubPointsCell = cellCub->getNumPoints();
00322 int numCubPointsSide = sideCub->getNumPoints();
00323
00324
00325
00326 FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim);
00327 FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim);
00328 FieldContainer<double> cub_weights_cell(numCubPointsCell);
00329 FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim);
00330 FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim);
00331 FieldContainer<double> jacobian_det_cell(1, numCubPointsCell);
00332 FieldContainer<double> weighted_measure_cell(1, numCubPointsCell);
00333
00334 FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell);
00335 FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00336 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00337 FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim);
00338 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00339 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00340 FieldContainer<double> fe_matrix(1, numFields, numFields);
00341
00342 FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell);
00343 FieldContainer<double> rhs_and_soln_vector(1, numFields);
00344
00345
00346 unsigned numSides = 4;
00347 FieldContainer<double> cub_points_side(numCubPointsSide, sideDim);
00348 FieldContainer<double> cub_weights_side(numCubPointsSide);
00349 FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim);
00350 FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim);
00351 FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim);
00352 FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide);
00353 FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide);
00354
00355 FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide);
00356 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00357 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00358 FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide);
00359 FieldContainer<double> neumann_fields_per_side(1, numFields);
00360
00361
00362 FieldContainer<double> value_of_basis_at_interp_points_ref(numFields, numInterpPoints);
00363 FieldContainer<double> transformed_value_of_basis_at_interp_points_ref(1, numFields, numInterpPoints);
00364 FieldContainer<double> interpolant(1, numInterpPoints);
00365
00366 FieldContainer<int> ipiv(numFields);
00367
00368
00369
00370
00371
00372
00373 cellCub->getCubature(cub_points_cell, cub_weights_cell);
00374
00375
00376 CellTools<double>::setJacobian(jacobian_cell, cub_points_cell, cell_nodes, cell);
00377 CellTools<double>::setJacobianInv(jacobian_inv_cell, jacobian_cell);
00378 CellTools<double>::setJacobianDet(jacobian_det_cell, jacobian_cell);
00379
00380
00381 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
00382
00384
00385
00386 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
00387
00388
00389 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
00390 value_of_basis_at_cub_points_cell);
00391
00392
00393 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
00394 weighted_measure_cell,
00395 transformed_value_of_basis_at_cub_points_cell);
00396
00397
00398 FunctionSpaceTools::integrate<double>(fe_matrix,
00399 transformed_value_of_basis_at_cub_points_cell,
00400 weighted_transformed_value_of_basis_at_cub_points_cell,
00401 COMP_BLAS);
00403
00405
00406
00407 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
00408
00409
00410 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
00411 jacobian_inv_cell,
00412 grad_of_basis_at_cub_points_cell);
00413
00414
00415 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
00416 weighted_measure_cell,
00417 transformed_grad_of_basis_at_cub_points_cell);
00418
00419
00420 FunctionSpaceTools::integrate<double>(fe_matrix,
00421 transformed_grad_of_basis_at_cub_points_cell,
00422 weighted_transformed_grad_of_basis_at_cub_points_cell,
00423 COMP_BLAS,
00424 true);
00426
00428
00429
00430 CellTools<double>::mapToPhysicalFrame(cub_points_cell_physical, cub_points_cell, cell_nodes, cell);
00431
00432
00433 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order, z_order);
00434
00435
00436 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
00437 rhs_at_cub_points_cell_physical,
00438 weighted_transformed_value_of_basis_at_cub_points_cell,
00439 COMP_BLAS);
00440
00441
00442 sideCub->getCubature(cub_points_side, cub_weights_side);
00443 for (unsigned i=0; i<numSides; i++) {
00444
00445 CellTools<double>::mapToReferenceSubcell(cub_points_side_refcell, cub_points_side, sideDim, (int)i, cell);
00446 CellTools<double>::setJacobian(jacobian_side_refcell, cub_points_side_refcell, cell_nodes, cell);
00447 CellTools<double>::setJacobianDet(jacobian_det_side_refcell, jacobian_side_refcell);
00448
00449
00450 FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_side_refcell,
00451 jacobian_side_refcell,
00452 cub_weights_side,
00453 i,
00454 cell);
00455
00456
00457 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
00458
00459 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
00460 value_of_basis_at_cub_points_side_refcell);
00461
00462
00463 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00464 weighted_measure_side_refcell,
00465 transformed_value_of_basis_at_cub_points_side_refcell);
00466
00467
00468
00469 CellTools<double>::mapToPhysicalFrame(cub_points_side_physical, cub_points_side_refcell, cell_nodes, cell);
00470
00471 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
00472 cell, (int)i, x_order, y_order, z_order);
00473
00474 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
00475 neumann_data_at_cub_points_side_physical,
00476 weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00477 COMP_BLAS);
00478
00479
00480 RealSpaceTools<double>::add(rhs_and_soln_vector, neumann_fields_per_side);;
00481 }
00483
00485
00486 int info = 0;
00487 Teuchos::LAPACK<int, double> solver;
00488 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
00490
00492
00493
00494 basis->getValues(value_of_basis_at_interp_points_ref, interp_points_ref, OPERATOR_VALUE);
00495
00496 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points_ref,
00497 value_of_basis_at_interp_points_ref);
00498 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points_ref);
00500
00501
00502
00503 RealSpaceTools<double>::subtract(interpolant, exact_solution);
00504
00505 *outStream << "\nRelative norm-2 error between exact solution polynomial of order ("
00506 << x_order << ", " << y_order << ", " << z_order
00507 << ") and finite element interpolant of order " << basis_order << ": "
00508 << RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00509 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) << "\n";
00510
00511 if (RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00512 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) > zero) {
00513 *outStream << "\n\nPatch test failed for solution polynomial order ("
00514 << x_order << ", " << y_order << ", " << z_order << ") and basis order " << basis_order << "\n\n";
00515 errorFlag++;
00516 }
00517 }
00518 }
00519 }
00520 }
00521
00522 }
00523
00524 catch (std::logic_error err) {
00525 *outStream << err.what() << "\n\n";
00526 errorFlag = -1000;
00527 };
00528
00529 if (errorFlag != 0)
00530 std::cout << "End Result: TEST FAILED\n";
00531 else
00532 std::cout << "End Result: TEST PASSED\n";
00533
00534
00535 std::cout.copyfmt(oldFormatState);
00536
00537 return errorFlag;
00538 }