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"
64 using namespace Intrepid;
70 const shards::CellTopology & ,
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);
127 const shards::CellTopology & parentCell,
128 int sideOrdinal,
int xd,
int yd,
int zd) {
130 int x = 0, y = 1, z = 2;
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);
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);
192 int 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_HEX_I2_FEM) |\n" \
215 <<
"| 1) Patch test involving mass and stiffness matrices, |\n" \
216 <<
"| for the Neumann problem on a physical parallelepiped |\n" \
217 <<
"| Omega with boundary Gamma. |\n" \
219 <<
"| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \
222 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
223 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
224 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
226 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
227 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
229 <<
"===============================================================================\n"\
230 <<
"| TEST 1: Patch test |\n"\
231 <<
"===============================================================================\n";
236 outStream -> precision(16);
243 shards::CellTopology cell(shards::getCellTopologyData< shards::Hexahedron<> >());
244 shards::CellTopology side(shards::getCellTopologyData< shards::Quadrilateral<> >());
245 int cellDim = cell.getDimension();
246 int sideDim = side.getDimension();
247 unsigned numSides = 6;
250 int numIntervals = 10;
251 int numInterpPoints = (numIntervals + 1)*(numIntervals + 1)*(numIntervals + 1);
254 for (
int k=0; k<=numIntervals; k++) {
255 for (
int j=0; j<=numIntervals; j++) {
256 for (
int i=0; i<=numIntervals; i++) {
257 interp_points_ref(counter,0) = i*(1.0/numIntervals)-1.0;
258 interp_points_ref(counter,1) = j*(1.0/numIntervals)-1.0;
259 interp_points_ref(counter,2) = k*(1.0/numIntervals)-1.0;
269 cell_nodes(0, 0, 0) = -5.0;
270 cell_nodes(0, 0, 1) = -1.0;
271 cell_nodes(0, 0, 2) = 0.0;
272 cell_nodes(0, 1, 0) = 4.0;
273 cell_nodes(0, 1, 1) = 1.0;
274 cell_nodes(0, 1, 2) = 1.0;
275 cell_nodes(0, 2, 0) = 8.0;
276 cell_nodes(0, 2, 1) = 3.0;
277 cell_nodes(0, 2, 2) = 1.0;
278 cell_nodes(0, 3, 0) = -1.0;
279 cell_nodes(0, 3, 1) = 1.0;
280 cell_nodes(0, 3, 2) = 0.0;
281 cell_nodes(0, 4, 0) = 5.0;
282 cell_nodes(0, 4, 1) = 9.0;
283 cell_nodes(0, 4, 2) = 1.0;
284 cell_nodes(0, 5, 0) = 14.0;
285 cell_nodes(0, 5, 1) = 11.0;
286 cell_nodes(0, 5, 2) = 2.0;
287 cell_nodes(0, 6, 0) = 18.0;
288 cell_nodes(0, 6, 1) = 13.0;
289 cell_nodes(0, 6, 2) = 2.0;
290 cell_nodes(0, 7, 0) = 9.0;
291 cell_nodes(0, 7, 1) = 11.0;
292 cell_nodes(0, 7, 2) = 1.0;
322 interp_points.resize(numInterpPoints, cellDim);
324 for (
int x_order=0; x_order <= max_order; x_order++) {
325 for (
int y_order=0; y_order <= max_order - x_order; y_order++) {
326 for (
int z_order=0; z_order <= max_order - x_order - y_order; z_order++) {
330 u_exact(exact_solution, interp_points, x_order, y_order, z_order);
335 double zero = basis_order*basis_order*basis_order*100*INTREPID_TOL;
338 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
340 int numFields = basis->getCardinality();
343 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.
create(cell, 2*basis_order);
344 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.
create(side, 2*basis_order);
345 int numCubPointsCell = cellCub->getNumPoints();
346 int numCubPointsSide = sideCub->getNumPoints();
360 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
362 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
363 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
379 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
380 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
396 cellCub->getCubature(cub_points_cell, cub_weights_cell);
404 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
409 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
412 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
413 value_of_basis_at_cub_points_cell);
416 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
417 weighted_measure_cell,
418 transformed_value_of_basis_at_cub_points_cell);
421 FunctionSpaceTools::integrate<double>(fe_matrix,
422 transformed_value_of_basis_at_cub_points_cell,
423 weighted_transformed_value_of_basis_at_cub_points_cell,
430 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
433 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
435 grad_of_basis_at_cub_points_cell);
438 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
439 weighted_measure_cell,
440 transformed_grad_of_basis_at_cub_points_cell);
443 FunctionSpaceTools::integrate<double>(fe_matrix,
444 transformed_grad_of_basis_at_cub_points_cell,
445 weighted_transformed_grad_of_basis_at_cub_points_cell,
456 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order, z_order);
459 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
460 rhs_at_cub_points_cell_physical,
461 weighted_transformed_value_of_basis_at_cub_points_cell,
465 sideCub->getCubature(cub_points_side, cub_weights_side);
466 for (
unsigned i=0; i<numSides; i++) {
473 FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_side_refcell,
474 jacobian_side_refcell,
480 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
482 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
483 value_of_basis_at_cub_points_side_refcell);
486 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
487 weighted_measure_side_refcell,
488 transformed_value_of_basis_at_cub_points_side_refcell);
494 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
495 cell, (
int)i, x_order, y_order, z_order);
497 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
498 neumann_data_at_cub_points_side_physical,
499 weighted_transformed_value_of_basis_at_cub_points_side_refcell,
510 Teuchos::LAPACK<int, double> solver;
511 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
517 basis->getValues(value_of_basis_at_interp_points_ref, interp_points_ref, OPERATOR_VALUE);
519 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points_ref,
520 value_of_basis_at_interp_points_ref);
521 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points_ref);
528 *outStream <<
"\nRelative norm-2 error between exact solution polynomial of order ("
529 << x_order <<
", " << y_order <<
", " << z_order
530 <<
") and finite element interpolant of order " << basis_order <<
": "
536 *outStream <<
"\n\nPatch test failed for solution polynomial order ("
537 << x_order <<
", " << y_order <<
", " << z_order <<
") and basis order " << basis_order <<
"\n\n";
546 catch (
const std::logic_error & err) {
547 *outStream << err.what() <<
"\n\n";
552 std::cout <<
"End Result: TEST FAILED\n";
554 std::cout <<
"End Result: TEST PASSED\n";
557 std::cout.copyfmt(oldFormatState);
Implementation of the serendipity-family H(grad)-compatible FEM basis of degree 2 on a Hexahedron cel...
int dimension(const int whichDim) const
Returns the specified dimension.
Header file for utility class to provide multidimensional containers.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Header file for the Intrepid::HGRAD_HEX_I2_FEM class.
A factory class that generates specific instances of cubatures.
Teuchos::RCP< Cubature< Scalar, ArrayPoint, ArrayWeight > > create(const shards::CellTopology &cellTopology, const std::vector< int > °ree)
Factory method.