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 & ,
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) * std::pow(points(cell,pt,y), yd);
94 for (
int cell=0; cell<result.
dimension(0); cell++) {
95 for (
int pt=0; pt<result.
dimension(1); pt++) {
96 result(cell,pt) -= yd*(yd-1)*std::pow(points(cell,pt,y), yd-2) * std::pow(points(cell,pt,x), xd);
102 for (
int cell=0; cell<result.
dimension(0); cell++) {
103 for (
int pt=0; pt<result.
dimension(1); pt++) {
104 result(cell,pt) += std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
115 const shards::CellTopology & parentCell,
116 int sideOrdinal,
int xd,
int yd) {
129 for (
int cell=0; cell<numCells; cell++) {
130 for (
int pt=0; pt<numPoints; pt++) {
131 grad_u(cell,pt,x) = xd*std::pow(points(cell,pt,x), xd-1) * std::pow(points(cell,pt,y), yd);
138 for (
int cell=0; cell<numCells; cell++) {
139 for (
int pt=0; pt<numPoints; pt++) {
140 grad_u(cell,pt,y) = yd*std::pow(points(cell,pt,y), yd-1) * std::pow(points(cell,pt,x), xd);
149 FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals,
true);
151 FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals);
158 for (
int cell=0; cell<result.
dimension(0); cell++) {
159 for (
int pt=0; pt<result.
dimension(1); pt++) {
160 result(cell,pt) = std::pow(points(pt,x), xd)*std::pow(points(pt,y), yd);
168 int main(
int argc,
char *argv[]) {
170 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
174 int iprint = argc - 1;
175 Teuchos::RCP<std::ostream> outStream;
176 Teuchos::oblackholestream bhs;
178 outStream = Teuchos::rcp(&std::cout,
false);
180 outStream = Teuchos::rcp(&bhs,
false);
183 Teuchos::oblackholestream oldFormatState;
184 oldFormatState.copyfmt(std::cout);
187 <<
"===============================================================================\n" \
189 <<
"| Unit Test (Basis_HGRAD_QUAD_Cn_FEM) |\n" \
191 <<
"| 1) Patch test involving mass and stiffness matrices, |\n" \
192 <<
"| for the Neumann problem on a physical parallelogram |\n" \
193 <<
"| AND a reference quad Omega with boundary Gamma. |\n" \
195 <<
"| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \
197 <<
"| For a generic parallelogram, the basis recovers a complete |\n" \
198 <<
"| polynomial space of order 2. On a (scaled and/or translated) |\n" \
199 <<
"| reference quad, the basis recovers a complete tensor product |\n" \
200 <<
"| space of order 2 (i.e. incl. the x^2*y, x*y^2, x^2*y^2 terms). |\n" \
202 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
203 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
204 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
206 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
207 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
209 <<
"===============================================================================\n"\
210 <<
"| TEST 1: Patch test |\n"\
211 <<
"===============================================================================\n";
216 outStream -> precision(16);
223 shards::CellTopology cell(shards::getCellTopologyData< shards::Quadrilateral<> >());
224 shards::CellTopology side(shards::getCellTopologyData< shards::Line<> >());
225 int cellDim = cell.getDimension();
226 int sideDim = side.getDimension();
229 int numIntervals = 10;
230 int numInterpPoints = (numIntervals + 1)*(numIntervals + 1);
233 for (
int j=0; j<=numIntervals; j++) {
234 for (
int i=0; i<=numIntervals; i++) {
235 interp_points_ref(counter,0) = i*(2.0/numIntervals)-1.0;
236 interp_points_ref(counter,1) = j*(2.0/numIntervals)-1.0;
243 cell_nodes[0].
resize(1, 4, cellDim);
244 cell_nodes[1].
resize(1, 4, cellDim);
247 cell_nodes[0](0, 0, 0) = -5.0;
248 cell_nodes[0](0, 0, 1) = -1.0;
249 cell_nodes[0](0, 1, 0) = 4.0;
250 cell_nodes[0](0, 1, 1) = 1.0;
251 cell_nodes[0](0, 2, 0) = 8.0;
252 cell_nodes[0](0, 2, 1) = 3.0;
253 cell_nodes[0](0, 3, 0) = -1.0;
254 cell_nodes[0](0, 3, 1) = 1.0;
256 cell_nodes[1](0, 0, 0) = -1.0;
257 cell_nodes[1](0, 0, 1) = -1.0;
258 cell_nodes[1](0, 1, 0) = 1.0;
259 cell_nodes[1](0, 1, 1) = -1.0;
260 cell_nodes[1](0, 2, 0) = 1.0;
261 cell_nodes[1](0, 2, 1) = 1.0;
262 cell_nodes[1](0, 3, 0) = -1.0;
263 cell_nodes[1](0, 3, 1) = 1.0;
265 std::stringstream mystream[2];
266 mystream[0].str(
"\n>> Now testing basis on a generic parallelogram ...\n");
267 mystream[1].str(
"\n>> Now testing basis on the reference quad ...\n");
269 for (
int pcell = 0; pcell < 2; pcell++) {
270 *outStream << mystream[pcell].str();
273 interp_points.resize(numInterpPoints, cellDim);
275 for (
int x_order=0; x_order <= max_order; x_order++) {
276 int max_y_order = max_order;
278 max_y_order -= x_order;
280 for (
int y_order=0; y_order <= max_y_order; y_order++) {
284 u_exact(exact_solution, interp_points, x_order, y_order);
289 double zero = basis_order*basis_order*100*INTREPID_TOL;
292 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
294 int numFields = basis->getCardinality();
297 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.
create(cell, 2*basis_order);
298 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.
create(side, 2*basis_order);
299 int numCubPointsCell = cellCub->getNumPoints();
300 int numCubPointsSide = sideCub->getNumPoints();
314 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
316 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
317 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
324 unsigned numSides = 4;
334 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
335 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
351 cellCub->getCubature(cub_points_cell, cub_weights_cell);
359 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
364 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
367 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
368 value_of_basis_at_cub_points_cell);
371 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
372 weighted_measure_cell,
373 transformed_value_of_basis_at_cub_points_cell);
376 FunctionSpaceTools::integrate<double>(fe_matrix,
377 transformed_value_of_basis_at_cub_points_cell,
378 weighted_transformed_value_of_basis_at_cub_points_cell,
385 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
388 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
390 grad_of_basis_at_cub_points_cell);
393 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
394 weighted_measure_cell,
395 transformed_grad_of_basis_at_cub_points_cell);
398 FunctionSpaceTools::integrate<double>(fe_matrix,
399 transformed_grad_of_basis_at_cub_points_cell,
400 weighted_transformed_grad_of_basis_at_cub_points_cell,
411 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order);
414 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
415 rhs_at_cub_points_cell_physical,
416 weighted_transformed_value_of_basis_at_cub_points_cell,
420 sideCub->getCubature(cub_points_side, cub_weights_side);
421 for (
unsigned i=0; i<numSides; i++) {
428 FunctionSpaceTools::computeEdgeMeasure<double>(weighted_measure_side_refcell,
429 jacobian_side_refcell,
435 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
437 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
438 value_of_basis_at_cub_points_side_refcell);
441 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
442 weighted_measure_side_refcell,
443 transformed_value_of_basis_at_cub_points_side_refcell);
449 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
450 cell, (
int)i, x_order, y_order);
452 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
453 neumann_data_at_cub_points_side_physical,
454 weighted_transformed_value_of_basis_at_cub_points_side_refcell,
465 Teuchos::LAPACK<int, double> solver;
466 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
472 basis->getValues(value_of_basis_at_interp_points, interp_points_ref, OPERATOR_VALUE);
474 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points,
475 value_of_basis_at_interp_points);
476 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points);
483 *outStream <<
"\nRelative norm-2 error between exact solution polynomial of order ("
484 << x_order <<
", " << y_order <<
") and finite element interpolant of order " << basis_order <<
": "
490 *outStream <<
"\n\nPatch test failed for solution polynomial order ("
491 << x_order <<
", " << y_order <<
") and basis order " << basis_order <<
"\n\n";
500 catch (
const std::logic_error & err) {
501 *outStream << err.what() <<
"\n\n";
506 std::cout <<
"End Result: TEST FAILED\n";
508 std::cout <<
"End Result: TEST PASSED\n";
511 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 abstract base class Intrepid::DefaultCubatureFactory.
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0...
Header file for the Intrepid::HGRAD_QUAD_Cn_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.