57 #include "Teuchos_oblackholestream.hpp"
58 #include "Teuchos_RCP.hpp"
59 #include "Teuchos_GlobalMPISession.hpp"
60 #include "Teuchos_SerialDenseMatrix.hpp"
61 #include "Teuchos_SerialDenseVector.hpp"
62 #include "Teuchos_LAPACK.hpp"
65 using namespace Intrepid;
71 const shards::CellTopology & ,
82 int x = 0, y = 1, z = 2;
86 for (
int cell=0; cell<result.
dimension(0); cell++) {
87 for (
int pt=0; pt<result.
dimension(1); pt++) {
88 result(cell,pt) = - xd*(xd-1)*std::pow(points(cell,pt,x), xd-2) *
89 std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
96 for (
int cell=0; cell<result.
dimension(0); cell++) {
97 for (
int pt=0; pt<result.
dimension(1); pt++) {
98 result(cell,pt) -= yd*(yd-1)*std::pow(points(cell,pt,y), yd-2) *
99 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd);
106 for (
int cell=0; cell<result.
dimension(0); cell++) {
107 for (
int pt=0; pt<result.
dimension(1); pt++) {
108 result(cell,pt) -= zd*(zd-1)*std::pow(points(cell,pt,z), zd-2) *
109 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
115 for (
int cell=0; cell<result.
dimension(0); cell++) {
116 for (
int pt=0; pt<result.
dimension(1); pt++) {
117 result(cell,pt) += std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
128 const shards::CellTopology & parentCell,
129 int sideOrdinal,
int xd,
int yd,
int zd) {
131 int x = 0, y = 1, z = 2;
142 for (
int cell=0; cell<numCells; cell++) {
143 for (
int pt=0; pt<numPoints; pt++) {
144 grad_u(cell,pt,x) = xd*std::pow(points(cell,pt,x), xd-1) *
145 std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
152 for (
int cell=0; cell<numCells; cell++) {
153 for (
int pt=0; pt<numPoints; pt++) {
154 grad_u(cell,pt,y) = yd*std::pow(points(cell,pt,y), yd-1) *
155 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd);
162 for (
int cell=0; cell<numCells; cell++) {
163 for (
int pt=0; pt<numPoints; pt++) {
164 grad_u(cell,pt,z) = zd*std::pow(points(cell,pt,z), zd-1) *
165 std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
174 FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals,
true);
176 FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals);
182 int x = 0, y = 1, z = 2;
183 for (
int cell=0; cell<result.
dimension(0); cell++) {
184 for (
int pt=0; pt<result.
dimension(1); pt++) {
185 result(cell,pt) = std::pow(points(pt,x), xd)*std::pow(points(pt,y), yd)*std::pow(points(pt,z), zd);
193 int main(
int argc,
char *argv[]) {
195 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
199 int iprint = argc - 1;
200 Teuchos::RCP<std::ostream> outStream;
201 Teuchos::oblackholestream bhs;
203 outStream = Teuchos::rcp(&std::cout,
false);
205 outStream = Teuchos::rcp(&bhs,
false);
208 Teuchos::oblackholestream oldFormatState;
209 oldFormatState.copyfmt(std::cout);
212 <<
"===============================================================================\n" \
214 <<
"| Unit Test (Basis_HGRAD_HEX_Cn_FEM) |\n" \
216 <<
"| 1) Patch test involving mass and stiffness matrices, |\n" \
217 <<
"| for the Neumann problem on a physical parallelepiped |\n" \
218 <<
"| AND a reference hex Omega with boundary Gamma. |\n" \
220 <<
"| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \
222 <<
"| For a generic parallelepiped, the basis recovers a complete |\n" \
223 <<
"| polynomial space of order 2. On a (scaled and/or translated) |\n" \
224 <<
"| reference hex, the basis recovers a complete tensor product |\n" \
225 <<
"| space of order 1 (i.e. incl. cross terms, e.g. x^2*y^2*z^2). |\n" \
227 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
228 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
229 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
231 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
232 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
234 <<
"===============================================================================\n"\
235 <<
"| TEST 1: Patch test |\n"\
236 <<
"===============================================================================\n";
241 outStream -> precision(16);
248 shards::CellTopology cell(shards::getCellTopologyData< shards::Hexahedron<> >());
249 shards::CellTopology side(shards::getCellTopologyData< shards::Quadrilateral<> >());
250 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());
251 int cellDim = cell.getDimension();
252 int sideDim = side.getDimension();
253 unsigned numSides = 6;
256 int numIntervals = 10;
257 int numInterpPoints = (numIntervals + 1)*(numIntervals + 1)*(numIntervals + 1);
260 for (
int k=0; k<=numIntervals; k++) {
261 for (
int j=0; j<=numIntervals; j++) {
262 for (
int i=0; i<=numIntervals; i++) {
263 interp_points_ref(counter,0) = i*(1.0/numIntervals)-1.0;
264 interp_points_ref(counter,1) = j*(1.0/numIntervals)-1.0;
265 interp_points_ref(counter,2) = k*(1.0/numIntervals)-1.0;
273 cell_nodes[0].
resize(1, 8, cellDim);
274 cell_nodes[1].
resize(1, 8, cellDim);
277 cell_nodes[0](0, 0, 0) = -5.0;
278 cell_nodes[0](0, 0, 1) = -1.0;
279 cell_nodes[0](0, 0, 2) = 0.0;
280 cell_nodes[0](0, 1, 0) = 4.0;
281 cell_nodes[0](0, 1, 1) = 1.0;
282 cell_nodes[0](0, 1, 2) = 1.0;
283 cell_nodes[0](0, 2, 0) = 8.0;
284 cell_nodes[0](0, 2, 1) = 3.0;
285 cell_nodes[0](0, 2, 2) = 1.0;
286 cell_nodes[0](0, 3, 0) = -1.0;
287 cell_nodes[0](0, 3, 1) = 1.0;
288 cell_nodes[0](0, 3, 2) = 0.0;
289 cell_nodes[0](0, 4, 0) = 5.0;
290 cell_nodes[0](0, 4, 1) = 9.0;
291 cell_nodes[0](0, 4, 2) = 1.0;
292 cell_nodes[0](0, 5, 0) = 14.0;
293 cell_nodes[0](0, 5, 1) = 11.0;
294 cell_nodes[0](0, 5, 2) = 2.0;
295 cell_nodes[0](0, 6, 0) = 18.0;
296 cell_nodes[0](0, 6, 1) = 13.0;
297 cell_nodes[0](0, 6, 2) = 2.0;
298 cell_nodes[0](0, 7, 0) = 9.0;
299 cell_nodes[0](0, 7, 1) = 11.0;
300 cell_nodes[0](0, 7, 2) = 1.0;
302 cell_nodes[1](0, 0, 0) = -1.0;
303 cell_nodes[1](0, 0, 1) = -1.0;
304 cell_nodes[1](0, 0, 2) = -1.0;
305 cell_nodes[1](0, 1, 0) = 1.0;
306 cell_nodes[1](0, 1, 1) = -1.0;
307 cell_nodes[1](0, 1, 2) = -1.0;
308 cell_nodes[1](0, 2, 0) = 1.0;
309 cell_nodes[1](0, 2, 1) = 1.0;
310 cell_nodes[1](0, 2, 2) = -1.0;
311 cell_nodes[1](0, 3, 0) = -1.0;
312 cell_nodes[1](0, 3, 1) = 1.0;
313 cell_nodes[1](0, 3, 2) = -1.0;
314 cell_nodes[1](0, 4, 0) = -1.0;
315 cell_nodes[1](0, 4, 1) = -1.0;
316 cell_nodes[1](0, 4, 2) = 1.0;
317 cell_nodes[1](0, 5, 0) = 1.0;
318 cell_nodes[1](0, 5, 1) = -1.0;
319 cell_nodes[1](0, 5, 2) = 1.0;
320 cell_nodes[1](0, 6, 0) = 1.0;
321 cell_nodes[1](0, 6, 1) = 1.0;
322 cell_nodes[1](0, 6, 2) = 1.0;
323 cell_nodes[1](0, 7, 0) = -1.0;
324 cell_nodes[1](0, 7, 1) = 1.0;
325 cell_nodes[1](0, 7, 2) = 1.0;
327 std::stringstream mystream[2];
328 mystream[0].str(
"\n>> Now testing basis on a generic parallelepiped ...\n");
329 mystream[1].str(
"\n>> Now testing basis on the reference hex ...\n");
332 for (
int pcell = 0; pcell < 2; pcell++) {
333 *outStream << mystream[pcell].str();
336 interp_points.resize(numInterpPoints, cellDim);
338 for (
int x_order=0; x_order <= max_order; x_order++) {
339 int max_y_order = max_order;
341 max_y_order -= x_order;
343 for (
int y_order=0; y_order <= max_y_order; y_order++) {
344 int max_z_order = max_order;
346 max_z_order -= x_order;
347 max_z_order -= y_order;
349 for (
int z_order=0; z_order <= max_z_order; z_order++) {
353 u_exact(exact_solution, interp_points, x_order, y_order, z_order);
355 int basis_order = max_order;
358 double zero = basis_order*basis_order*basis_order*100*INTREPID_TOL;
362 PointTools::getLattice<double,FieldContainer<double> >(pts,line,basis_order);
364 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
366 int numFields = basis->getCardinality();
369 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.
create(cell, 2*basis_order);
370 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.
create(side, 2*basis_order);
371 int numCubPointsCell = cellCub->getNumPoints();
372 int numCubPointsSide = sideCub->getNumPoints();
386 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
388 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
389 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
405 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
406 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
422 cellCub->getCubature(cub_points_cell, cub_weights_cell);
430 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
435 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
438 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
439 value_of_basis_at_cub_points_cell);
442 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
443 weighted_measure_cell,
444 transformed_value_of_basis_at_cub_points_cell);
447 FunctionSpaceTools::integrate<double>(fe_matrix,
448 transformed_value_of_basis_at_cub_points_cell,
449 weighted_transformed_value_of_basis_at_cub_points_cell,
456 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
459 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
461 grad_of_basis_at_cub_points_cell);
464 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
465 weighted_measure_cell,
466 transformed_grad_of_basis_at_cub_points_cell);
469 FunctionSpaceTools::integrate<double>(fe_matrix,
470 transformed_grad_of_basis_at_cub_points_cell,
471 weighted_transformed_grad_of_basis_at_cub_points_cell,
482 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order, z_order);
485 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
486 rhs_at_cub_points_cell_physical,
487 weighted_transformed_value_of_basis_at_cub_points_cell,
491 sideCub->getCubature(cub_points_side, cub_weights_side);
492 for (
unsigned i=0; i<numSides; i++) {
499 FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_side_refcell,
500 jacobian_side_refcell,
506 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
508 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
509 value_of_basis_at_cub_points_side_refcell);
512 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
513 weighted_measure_side_refcell,
514 transformed_value_of_basis_at_cub_points_side_refcell);
520 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
521 cell, (
int)i, x_order, y_order, z_order);
523 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
524 neumann_data_at_cub_points_side_physical,
525 weighted_transformed_value_of_basis_at_cub_points_side_refcell,
536 Teuchos::LAPACK<int, double> solver;
537 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
543 basis->getValues(value_of_basis_at_interp_points_ref, interp_points_ref, OPERATOR_VALUE);
545 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points_ref,
546 value_of_basis_at_interp_points_ref);
547 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points_ref);
554 *outStream <<
"\nRelative norm-2 error between exact solution polynomial of order ("
555 << x_order <<
", " << y_order <<
", " << z_order
556 <<
") and finite element interpolant of order " << basis_order <<
": "
562 *outStream <<
"\n\nPatch test failed for solution polynomial order ("
563 << x_order <<
", " << y_order <<
", " << z_order <<
") and basis order " << basis_order <<
"\n\n";
573 catch (
const std::logic_error & err) {
574 *outStream << err.what() <<
"\n\n";
579 std::cout <<
"End Result: TEST FAILED\n";
581 std::cout <<
"End Result: TEST PASSED\n";
584 std::cout.copyfmt(oldFormatState);
int dimension(const int whichDim) const
Returns the specified dimension.
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HGRAD_HEX_Cn_FEM class.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0...
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.