50#include "Intrepid_HGRAD_TET_C1_FEM.hpp"
56#include "Teuchos_oblackholestream.hpp"
57#include "Teuchos_RCP.hpp"
58#include "Teuchos_GlobalMPISession.hpp"
59#include "Teuchos_SerialDenseMatrix.hpp"
60#include "Teuchos_SerialDenseVector.hpp"
61#include "Teuchos_LAPACK.hpp"
64using namespace Intrepid;
66void rhsFunc(FieldContainer<double> &,
const FieldContainer<double> &,
int,
int,
int);
67void neumann(FieldContainer<double> & ,
68 const FieldContainer<double> & ,
69 const FieldContainer<double> & ,
70 const shards::CellTopology & ,
72void u_exact(FieldContainer<double> &,
const FieldContainer<double> &,
int,
int,
int);
75void rhsFunc(FieldContainer<double> & result,
76 const FieldContainer<double> & points,
81 int x = 0, y = 1, z = 2;
85 for (
int cell=0; cell<result.dimension(0); cell++) {
86 for (
int pt=0; pt<result.dimension(1); pt++) {
87 result(cell,pt) = - xd*(xd-1)*std::pow(points(cell,pt,x), xd-2) *
88 std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
95 for (
int cell=0; cell<result.dimension(0); cell++) {
96 for (
int pt=0; pt<result.dimension(1); pt++) {
97 result(cell,pt) -= yd*(yd-1)*std::pow(points(cell,pt,y), yd-2) *
98 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd);
105 for (
int cell=0; cell<result.dimension(0); cell++) {
106 for (
int pt=0; pt<result.dimension(1); pt++) {
107 result(cell,pt) -= zd*(zd-1)*std::pow(points(cell,pt,z), zd-2) *
108 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
114 for (
int cell=0; cell<result.dimension(0); cell++) {
115 for (
int pt=0; pt<result.dimension(1); pt++) {
116 result(cell,pt) += std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
125 const FieldContainer<double> & points,
126 const FieldContainer<double> & jacs,
127 const shards::CellTopology & parentCell,
128 int sideOrdinal,
int xd,
int yd,
int zd) {
130 int x = 0, y = 1, z = 2;
132 int numCells = result.dimension(0);
133 int numPoints = result.dimension(1);
135 FieldContainer<double> grad_u(numCells, numPoints, 3);
136 FieldContainer<double> side_normals(numCells, numPoints, 3);
137 FieldContainer<double> normal_lengths(numCells, numPoints);
141 for (
int cell=0; cell<numCells; cell++) {
142 for (
int pt=0; pt<numPoints; pt++) {
143 grad_u(cell,pt,x) = xd*std::pow(points(cell,pt,x), xd-1) *
144 std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
151 for (
int cell=0; cell<numCells; cell++) {
152 for (
int pt=0; pt<numPoints; pt++) {
153 grad_u(cell,pt,y) = yd*std::pow(points(cell,pt,y), yd-1) *
154 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd);
161 for (
int cell=0; cell<numCells; cell++) {
162 for (
int pt=0; pt<numPoints; pt++) {
163 grad_u(cell,pt,z) = zd*std::pow(points(cell,pt,z), zd-1) *
164 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
173 FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals,
true);
175 FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals);
180void u_exact(FieldContainer<double> & result,
const FieldContainer<double> & points,
int xd,
int yd,
int zd) {
181 int x = 0, y = 1, z = 2;
182 for (
int cell=0; cell<result.dimension(0); cell++) {
183 for (
int pt=0; pt<result.dimension(1); pt++) {
184 result(cell,pt) = std::pow(points(pt,x), xd)*std::pow(points(pt,y), yd)*std::pow(points(pt,z), zd);
192int main(
int argc,
char *argv[]) {
194 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
198 int iprint = argc - 1;
199 Teuchos::RCP<std::ostream> outStream;
200 Teuchos::oblackholestream bhs;
202 outStream = Teuchos::rcp(&std::cout,
false);
204 outStream = Teuchos::rcp(&bhs,
false);
207 Teuchos::oblackholestream oldFormatState;
208 oldFormatState.copyfmt(std::cout);
211 <<
"===============================================================================\n" \
213 <<
"| Unit Test (Basis_HGRAD_TET_C1_FEM) |\n" \
215 <<
"| 1) Patch test involving mass and stiffness matrices, |\n" \
216 <<
"| for the Neumann problem on a tetrahedral patch |\n" \
217 <<
"| Omega with boundary Gamma. |\n" \
219 <<
"| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \
221 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
222 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
223 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
225 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
226 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
228 <<
"===============================================================================\n"\
229 <<
"| TEST 1: Patch test |\n"\
230 <<
"===============================================================================\n";
235 outStream -> precision(16);
241 DefaultCubatureFactory<double> cubFactory;
242 shards::CellTopology cell(shards::getCellTopologyData< shards::Tetrahedron<> >());
243 shards::CellTopology side(shards::getCellTopologyData< shards::Triangle<> >());
244 int cellDim = cell.getDimension();
245 int sideDim = side.getDimension();
248 int numIntervals = 10;
249 int numInterpPoints = ((numIntervals + 1)*(numIntervals + 2)*(numIntervals + 3))/6;
250 FieldContainer<double> interp_points_ref(numInterpPoints, 3);
252 for (
int k=0; k<=numIntervals; k++) {
253 for (
int j=0; j<=numIntervals; j++) {
254 for (
int i=0; i<=numIntervals; i++) {
255 if (i+j+k <= numIntervals) {
256 interp_points_ref(counter,0) = i*(1.0/numIntervals);
257 interp_points_ref(counter,1) = j*(1.0/numIntervals);
258 interp_points_ref(counter,2) = k*(1.0/numIntervals);
266 FieldContainer<double> cell_nodes(1, 4, cellDim);
268 cell_nodes(0, 0, 0) = -1.0;
269 cell_nodes(0, 0, 1) = -2.0;
270 cell_nodes(0, 0, 2) = 0.0;
271 cell_nodes(0, 1, 0) = 6.0;
272 cell_nodes(0, 1, 1) = 2.0;
273 cell_nodes(0, 1, 2) = 0.0;
274 cell_nodes(0, 2, 0) = -5.0;
275 cell_nodes(0, 2, 1) = 1.0;
276 cell_nodes(0, 2, 2) = 0.0;
277 cell_nodes(0, 3, 0) = -4.0;
278 cell_nodes(0, 3, 1) = -1.0;
279 cell_nodes(0, 3, 2) = 3.0;
307 FieldContainer<double> interp_points(1, numInterpPoints, cellDim);
309 interp_points.resize(numInterpPoints, cellDim);
311 for (
int x_order=0; x_order <= max_order; x_order++) {
312 for (
int y_order=0; y_order <= max_order-x_order; y_order++) {
313 for (
int z_order=0; z_order <= max_order-x_order-y_order; z_order++) {
316 FieldContainer<double> exact_solution(1, numInterpPoints);
317 u_exact(exact_solution, interp_points, x_order, y_order, z_order);
322 double zero = basis_order*basis_order*basis_order*100*INTREPID_TOL;
325 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
327 int numFields = basis->getCardinality();
330 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*basis_order);
331 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*basis_order);
332 int numCubPointsCell = cellCub->getNumPoints();
333 int numCubPointsSide = sideCub->getNumPoints();
337 FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim);
338 FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim);
339 FieldContainer<double> cub_weights_cell(numCubPointsCell);
340 FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim);
341 FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim);
342 FieldContainer<double> jacobian_det_cell(1, numCubPointsCell);
343 FieldContainer<double> weighted_measure_cell(1, numCubPointsCell);
345 FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell);
346 FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
347 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
348 FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim);
349 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
350 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
351 FieldContainer<double> fe_matrix(1, numFields, numFields);
353 FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell);
354 FieldContainer<double> rhs_and_soln_vector(1, numFields);
357 unsigned numSides = 4;
358 FieldContainer<double> cub_points_side(numCubPointsSide, sideDim);
359 FieldContainer<double> cub_weights_side(numCubPointsSide);
360 FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim);
361 FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim);
362 FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim);
363 FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide);
364 FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide);
366 FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide);
367 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
368 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
369 FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide);
370 FieldContainer<double> neumann_fields_per_side(1, numFields);
373 FieldContainer<double> value_of_basis_at_interp_points_ref(numFields, numInterpPoints);
374 FieldContainer<double> transformed_value_of_basis_at_interp_points_ref(1, numFields, numInterpPoints);
375 FieldContainer<double> interpolant(1, numInterpPoints);
377 FieldContainer<int> ipiv(numFields);
384 cellCub->getCubature(cub_points_cell, cub_weights_cell);
392 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
397 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
400 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
401 value_of_basis_at_cub_points_cell);
404 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
405 weighted_measure_cell,
406 transformed_value_of_basis_at_cub_points_cell);
409 FunctionSpaceTools::integrate<double>(fe_matrix,
410 transformed_value_of_basis_at_cub_points_cell,
411 weighted_transformed_value_of_basis_at_cub_points_cell,
418 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
421 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
423 grad_of_basis_at_cub_points_cell);
426 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
427 weighted_measure_cell,
428 transformed_grad_of_basis_at_cub_points_cell);
431 FunctionSpaceTools::integrate<double>(fe_matrix,
432 transformed_grad_of_basis_at_cub_points_cell,
433 weighted_transformed_grad_of_basis_at_cub_points_cell,
444 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order, z_order);
447 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
448 rhs_at_cub_points_cell_physical,
449 weighted_transformed_value_of_basis_at_cub_points_cell,
453 sideCub->getCubature(cub_points_side, cub_weights_side);
454 for (
unsigned i=0; i<numSides; i++) {
461 FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_side_refcell,
462 jacobian_side_refcell,
468 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
470 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
471 value_of_basis_at_cub_points_side_refcell);
474 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
475 weighted_measure_side_refcell,
476 transformed_value_of_basis_at_cub_points_side_refcell);
482 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
483 cell, (
int)i, x_order, y_order, z_order);
485 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
486 neumann_data_at_cub_points_side_physical,
487 weighted_transformed_value_of_basis_at_cub_points_side_refcell,
498 Teuchos::LAPACK<int, double> solver;
499 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
505 basis->getValues(value_of_basis_at_interp_points_ref, interp_points_ref, OPERATOR_VALUE);
507 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points_ref,
508 value_of_basis_at_interp_points_ref);
509 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points_ref);
516 *outStream <<
"\nRelative norm-2 error between exact solution polynomial of order ("
517 << x_order <<
", " << y_order <<
", " << z_order
518 <<
") and finite element interpolant of order " << basis_order <<
": "
524 *outStream <<
"\n\nPatch test failed for solution polynomial order ("
525 << x_order <<
", " << y_order <<
", " << z_order <<
") and basis order " << basis_order <<
"\n\n";
534 catch (
const std::logic_error & err) {
535 *outStream << err.what() <<
"\n\n";
540 std::cout <<
"End Result: TEST FAILED\n";
542 std::cout <<
"End Result: TEST PASSED\n";
545 std::cout.copyfmt(oldFormatState);
void rhsFunc(FieldContainer< double > &, const FieldContainer< double > &, int, int, int)
right-hand side function
void u_exact(FieldContainer< double > &, const FieldContainer< double > &, int, int, int)
exact solution
void neumann(FieldContainer< double > &, const FieldContainer< double > &, const FieldContainer< double > &, const shards::CellTopology &, int, int, int, int)
neumann boundary conditions
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Header file for utility class to provide multidimensional containers.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell.