50 #include "Intrepid_HGRAD_WEDGE_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"
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_WEDGE_C1_FEM) |\n" \
215 <<
"| 1) Patch test involving mass and stiffness matrices, |\n" \
216 <<
"| for the Neumann problem on a wedge 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" \
224 <<
"| Mauro Perego (mperego@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::Wedge<> >());
244 shards::CellTopology sideQ(shards::getCellTopologyData< shards::Quadrilateral<> >());
245 shards::CellTopology sideT(shards::getCellTopologyData< shards::Triangle<> >());
246 int cellDim = cell.getDimension();
247 int sideQDim = sideQ.getDimension();
248 int sideTDim = sideT.getDimension();
251 int numIntervals = 10;
252 int numInterpPoints = ((numIntervals + 1)*(numIntervals + 2)*(numIntervals + 3))/6;
255 for (
int k=0; k<=numIntervals; k++) {
256 for (
int j=0; j<=numIntervals; j++) {
257 for (
int i=0; i<=numIntervals; i++) {
258 if (i+j+k <= numIntervals) {
259 interp_points_ref(counter,0) = i*(1.0/numIntervals);
260 interp_points_ref(counter,1) = j*(1.0/numIntervals);
261 interp_points_ref(counter,2) = k*(1.0/numIntervals);
271 cell_nodes(0, 0, 0) = -3.0;
272 cell_nodes(0, 0, 1) = -5.0;
273 cell_nodes(0, 0, 2) = -5.0;
274 cell_nodes(0, 1, 0) = -1.0;
275 cell_nodes(0, 1, 1) = -4.0;
276 cell_nodes(0, 1, 2) = -4.0;
277 cell_nodes(0, 2, 0) = 0.0;
278 cell_nodes(0, 2, 1) = 0.0;
279 cell_nodes(0, 2, 2) = -5.0;
280 cell_nodes(0, 3, 0) = 5.0;
281 cell_nodes(0, 3, 1) = -1.0;
282 cell_nodes(0, 3, 2) = 1.0;
283 cell_nodes(0, 4, 0) = 7.0;
284 cell_nodes(0, 4, 1) = 0.0;
285 cell_nodes(0, 4, 2) = 2.0;
286 cell_nodes(0, 5, 0) = 8.0;
287 cell_nodes(0, 5, 1) = 4.0;
288 cell_nodes(0, 5, 2) = 1.0;
333 interp_points.resize(numInterpPoints, cellDim);
335 for (
int x_order=0; x_order <= max_order; x_order++) {
336 for (
int y_order=0; y_order <= max_order-x_order; y_order++) {
337 for (
int z_order=0; z_order <= max_order-x_order-y_order; z_order++) {
341 u_exact(exact_solution, interp_points, x_order, y_order, z_order);
346 double zero = basis_order*basis_order*basis_order*100*INTREPID_TOL;
349 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
351 int numFields = basis->getCardinality();
354 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.
create(cell, 2*basis_order);
355 Teuchos::RCP<Cubature<double> > sideQCub = cubFactory.
create(sideQ, 2*basis_order);
356 Teuchos::RCP<Cubature<double> > sideTCub = cubFactory.
create(sideT, 2*basis_order);
357 int numCubPointsCell = cellCub->getNumPoints();
358 int numCubPointsSideQ = sideQCub->getNumPoints();
359 int numCubPointsSideT = sideTCub->getNumPoints();
373 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
375 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
376 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
383 unsigned numSides = 5;
384 unsigned numSidesQ = 3;
402 FieldContainer<double> transformed_value_of_basis_at_cub_points_sideQ_refcell(1, numFields, numCubPointsSideQ);
403 FieldContainer<double> transformed_value_of_basis_at_cub_points_sideT_refcell(1, numFields, numCubPointsSideT);
404 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_sideQ_refcell(1, numFields, numCubPointsSideQ);
405 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_sideT_refcell(1, numFields, numCubPointsSideT);
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 sideQCub->getCubature(cub_points_sideQ, cub_weights_sideQ);
492 sideTCub->getCubature(cub_points_sideT, cub_weights_sideT);
494 for (
unsigned i=0; i<numSidesQ; i++) {
501 FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_sideQ_refcell,
502 jacobian_sideQ_refcell,
508 basis->getValues(value_of_basis_at_cub_points_sideQ_refcell, cub_points_sideQ_refcell, OPERATOR_VALUE);
510 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_sideQ_refcell,
511 value_of_basis_at_cub_points_sideQ_refcell);
514 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_sideQ_refcell,
515 weighted_measure_sideQ_refcell,
516 transformed_value_of_basis_at_cub_points_sideQ_refcell);
522 neumann(neumann_data_at_cub_points_sideQ_physical, cub_points_sideQ_physical, jacobian_sideQ_refcell,
523 cell, (
int)i, x_order, y_order, z_order);
525 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
526 neumann_data_at_cub_points_sideQ_physical,
527 weighted_transformed_value_of_basis_at_cub_points_sideQ_refcell,
534 for (
unsigned i=numSidesQ; i<numSides; i++) {
541 FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_sideT_refcell,
542 jacobian_sideT_refcell,
548 basis->getValues(value_of_basis_at_cub_points_sideT_refcell, cub_points_sideT_refcell, OPERATOR_VALUE);
550 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_sideT_refcell,
551 value_of_basis_at_cub_points_sideT_refcell);
554 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_sideT_refcell,
555 weighted_measure_sideT_refcell,
556 transformed_value_of_basis_at_cub_points_sideT_refcell);
562 neumann(neumann_data_at_cub_points_sideT_physical, cub_points_sideT_physical, jacobian_sideT_refcell,
563 cell, (
int)i, x_order, y_order, z_order);
565 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
566 neumann_data_at_cub_points_sideT_physical,
567 weighted_transformed_value_of_basis_at_cub_points_sideT_refcell,
578 Teuchos::LAPACK<int, double> solver;
579 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
585 basis->getValues(value_of_basis_at_interp_points_ref, interp_points_ref, OPERATOR_VALUE);
587 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points_ref,
588 value_of_basis_at_interp_points_ref);
589 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points_ref);
596 *outStream <<
"\nRelative norm-2 error between exact solution polynomial of order ("
597 << x_order <<
", " << y_order <<
", " << z_order
598 <<
") and finite element interpolant of order " << basis_order <<
": "
604 *outStream <<
"\n\nPatch test failed for solution polynomial order ("
605 << x_order <<
", " << y_order <<
", " << z_order <<
") and basis order " << basis_order <<
"\n\n";
614 catch (
const std::logic_error & err) {
615 *outStream << err.what() <<
"\n\n";
620 std::cout <<
"End Result: TEST FAILED\n";
622 std::cout <<
"End Result: TEST PASSED\n";
625 std::cout.copyfmt(oldFormatState);
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
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.
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.