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);
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);
75void rhsFunc(FieldContainer<double> & result,
76 const FieldContainer<double> & points,
84 for (
int cell=0; cell<result.dimension(0); cell++) {
85 for (
int pt=0; pt<result.dimension(1); pt++) {
86 result(cell,pt) = - xd*(xd-1)*std::pow(points(cell,pt,x), xd-2) * std::pow(points(cell,pt,y), yd);
93 for (
int cell=0; cell<result.dimension(0); cell++) {
94 for (
int pt=0; pt<result.dimension(1); pt++) {
95 result(cell,pt) -= yd*(yd-1)*std::pow(points(cell,pt,y), yd-2) * std::pow(points(cell,pt,x), xd);
101 for (
int cell=0; cell<result.dimension(0); cell++) {
102 for (
int pt=0; pt<result.dimension(1); pt++) {
103 result(cell,pt) += std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
112 const FieldContainer<double> & points,
113 const FieldContainer<double> & jacs,
114 const shards::CellTopology & parentCell,
115 int sideOrdinal,
int xd,
int yd) {
119 int numCells = result.dimension(0);
120 int numPoints = result.dimension(1);
122 FieldContainer<double> grad_u(numCells, numPoints, 2);
123 FieldContainer<double> side_normals(numCells, numPoints, 2);
124 FieldContainer<double> normal_lengths(numCells, numPoints);
128 for (
int cell=0; cell<numCells; cell++) {
129 for (
int pt=0; pt<numPoints; pt++) {
130 grad_u(cell,pt,x) = xd*std::pow(points(cell,pt,x), xd-1) * std::pow(points(cell,pt,y), yd);
137 for (
int cell=0; cell<numCells; cell++) {
138 for (
int pt=0; pt<numPoints; pt++) {
139 grad_u(cell,pt,y) = yd*std::pow(points(cell,pt,y), yd-1) * std::pow(points(cell,pt,x), xd);
148 FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals,
true);
150 FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals);
155void u_exact(FieldContainer<double> & result,
const FieldContainer<double> & points,
int xd,
int yd) {
157 for (
int cell=0; cell<result.dimension(0); cell++) {
158 for (
int pt=0; pt<result.dimension(1); pt++) {
159 result(cell,pt) = std::pow(points(pt,x), xd)*std::pow(points(pt,y), yd);
167int main(
int argc,
char *argv[]) {
169 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
173 int iprint = argc - 1;
174 Teuchos::RCP<std::ostream> outStream;
175 Teuchos::oblackholestream bhs;
177 outStream = Teuchos::rcp(&std::cout,
false);
179 outStream = Teuchos::rcp(&bhs,
false);
182 Teuchos::oblackholestream oldFormatState;
183 oldFormatState.copyfmt(std::cout);
186 <<
"===============================================================================\n" \
188 <<
"| Unit Test (Basis_HGRAD_TRI_C1_FEM) |\n" \
190 <<
"| 1) Patch test involving mass and stiffness matrices, |\n" \
191 <<
"| for the Neumann problem on a triangular patch |\n" \
192 <<
"| Omega with boundary Gamma. |\n" \
194 <<
"| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \
196 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
197 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
198 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
200 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
201 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
203 <<
"===============================================================================\n"\
204 <<
"| TEST 1: Patch test |\n"\
205 <<
"===============================================================================\n";
210 outStream -> precision(16);
216 DefaultCubatureFactory<double> cubFactory;
217 shards::CellTopology cell(shards::getCellTopologyData< shards::Triangle<> >());
218 shards::CellTopology side(shards::getCellTopologyData< shards::Line<> >());
219 int cellDim = cell.getDimension();
220 int sideDim = side.getDimension();
223 int numIntervals = 10;
224 int numInterpPoints = ((numIntervals + 1)*(numIntervals + 2))/2;
225 FieldContainer<double> interp_points_ref(numInterpPoints, 2);
227 for (
int j=0; j<=numIntervals; j++) {
228 for (
int i=0; i<=numIntervals; i++) {
229 if (i <= numIntervals-j) {
230 interp_points_ref(counter,0) = i*(1.0/numIntervals);
231 interp_points_ref(counter,1) = j*(1.0/numIntervals);
238 FieldContainer<double> cell_nodes(1, 3, cellDim);
240 cell_nodes(0, 0, 0) = 0.1;
241 cell_nodes(0, 0, 1) = -0.1;
242 cell_nodes(0, 1, 0) = 1.1;
243 cell_nodes(0, 1, 1) = -0.1;
244 cell_nodes(0, 2, 0) = 0.1;
245 cell_nodes(0, 2, 1) = 0.9;
254 FieldContainer<double> interp_points(1, numInterpPoints, cellDim);
256 interp_points.resize(numInterpPoints, cellDim);
258 for (
int x_order=0; x_order <= max_order; x_order++) {
259 for (
int y_order=0; y_order <= max_order-x_order; y_order++) {
262 FieldContainer<double> exact_solution(1, numInterpPoints);
263 u_exact(exact_solution, interp_points, x_order, y_order);
268 double zero = basis_order*basis_order*100*INTREPID_TOL;
271 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
273 int numFields = basis->getCardinality();
276 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*basis_order);
277 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*basis_order);
278 int numCubPointsCell = cellCub->getNumPoints();
279 int numCubPointsSide = sideCub->getNumPoints();
283 FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim);
284 FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim);
285 FieldContainer<double> cub_weights_cell(numCubPointsCell);
286 FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim);
287 FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim);
288 FieldContainer<double> jacobian_det_cell(1, numCubPointsCell);
289 FieldContainer<double> weighted_measure_cell(1, numCubPointsCell);
291 FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell);
292 FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
293 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
294 FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim);
295 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
296 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
297 FieldContainer<double> fe_matrix(1, numFields, numFields);
299 FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell);
300 FieldContainer<double> rhs_and_soln_vector(1, numFields);
303 unsigned numSides = 3;
304 FieldContainer<double> cub_points_side(numCubPointsSide, sideDim);
305 FieldContainer<double> cub_weights_side(numCubPointsSide);
306 FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim);
307 FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim);
308 FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim);
309 FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide);
310 FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide);
312 FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide);
313 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
314 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
315 FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide);
316 FieldContainer<double> neumann_fields_per_side(1, numFields);
319 FieldContainer<double> value_of_basis_at_interp_points(numFields, numInterpPoints);
320 FieldContainer<double> transformed_value_of_basis_at_interp_points(1, numFields, numInterpPoints);
321 FieldContainer<double> interpolant(1, numInterpPoints);
323 FieldContainer<int> ipiv(numFields);
330 cellCub->getCubature(cub_points_cell, cub_weights_cell);
338 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
343 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
346 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
347 value_of_basis_at_cub_points_cell);
350 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
351 weighted_measure_cell,
352 transformed_value_of_basis_at_cub_points_cell);
355 FunctionSpaceTools::integrate<double>(fe_matrix,
356 transformed_value_of_basis_at_cub_points_cell,
357 weighted_transformed_value_of_basis_at_cub_points_cell,
364 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
367 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
369 grad_of_basis_at_cub_points_cell);
372 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
373 weighted_measure_cell,
374 transformed_grad_of_basis_at_cub_points_cell);
377 FunctionSpaceTools::integrate<double>(fe_matrix,
378 transformed_grad_of_basis_at_cub_points_cell,
379 weighted_transformed_grad_of_basis_at_cub_points_cell,
390 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order);
393 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
394 rhs_at_cub_points_cell_physical,
395 weighted_transformed_value_of_basis_at_cub_points_cell,
399 sideCub->getCubature(cub_points_side, cub_weights_side);
400 for (
unsigned i=0; i<numSides; i++) {
407 FunctionSpaceTools::computeEdgeMeasure<double>(weighted_measure_side_refcell,
408 jacobian_side_refcell,
414 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
416 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
417 value_of_basis_at_cub_points_side_refcell);
420 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
421 weighted_measure_side_refcell,
422 transformed_value_of_basis_at_cub_points_side_refcell);
428 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
429 cell, (
int)i, x_order, y_order);
431 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
432 neumann_data_at_cub_points_side_physical,
433 weighted_transformed_value_of_basis_at_cub_points_side_refcell,
444 Teuchos::LAPACK<int, double> solver;
445 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
451 basis->getValues(value_of_basis_at_interp_points, interp_points_ref, OPERATOR_VALUE);
453 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points,
454 value_of_basis_at_interp_points);
455 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points);
462 *outStream <<
"\nRelative norm-2 error between exact solution polynomial of order ("
463 << x_order <<
", " << y_order <<
") and finite element interpolant of order " << basis_order <<
": "
469 *outStream <<
"\n\nPatch test failed for solution polynomial order ("
470 << x_order <<
", " << y_order <<
") and basis order " << basis_order <<
"\n\n";
478 catch (
const std::logic_error & err) {
479 *outStream << err.what() <<
"\n\n";
484 std::cout <<
"End Result: TEST FAILED\n";
486 std::cout <<
"End Result: TEST PASSED\n";
489 std::cout.copyfmt(oldFormatState);
void neumann(FieldContainer< double > &, const FieldContainer< double > &, const FieldContainer< double > &, const shards::CellTopology &, int, int, int)
neumann boundary conditions
void u_exact(FieldContainer< double > &, const FieldContainer< double > &, int, int)
exact solution
void rhsFunc(FieldContainer< double > &, const FieldContainer< double > &, int, int)
right-hand side function
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HGRAD_TRI_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell.