Intrepid
test_02.cpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Pavel Bochev (pbboche@sandia.gov)
38// Denis Ridzal (dridzal@sandia.gov), or
39// Kara Peterson (kjpeter@sandia.gov)
40//
41// ************************************************************************
42// @HEADER
43
48
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"
63
64using namespace std;
65using namespace Intrepid;
66
71 const shards::CellTopology & ,
72 int, int, int);
74
77 const FieldContainer<double> & points,
78 int xd,
79 int yd) {
80
81 int x = 0, y = 1;
82
83 // second x-derivatives of u
84 if (xd > 1) {
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);
88 }
89 }
90 }
91
92 // second y-derivatives of u
93 if (yd > 1) {
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);
97 }
98 }
99 }
100
101 // add u
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);
105 }
106 }
107
108}
109
110
113 const FieldContainer<double> & points,
114 const FieldContainer<double> & jacs,
115 const shards::CellTopology & parentCell,
116 int sideOrdinal, int xd, int yd) {
117
118 int x = 0, y = 1;
119
120 int numCells = result.dimension(0);
121 int numPoints = result.dimension(1);
122
123 FieldContainer<double> grad_u(numCells, numPoints, 2);
124 FieldContainer<double> side_normals(numCells, numPoints, 2);
125 FieldContainer<double> normal_lengths(numCells, numPoints);
126
127 // first x-derivatives of u
128 if (xd > 0) {
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);
132 }
133 }
134 }
135
136 // first y-derivatives of u
137 if (yd > 0) {
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);
141 }
142 }
143 }
144
145 CellTools<double>::getPhysicalSideNormals(side_normals, jacs, sideOrdinal, parentCell);
146
147 // scale normals
148 RealSpaceTools<double>::vectorNorm(normal_lengths, side_normals, NORM_TWO);
149 FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals, true);
150
151 FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals);
152
153}
154
156void u_exact(FieldContainer<double> & result, const FieldContainer<double> & points, int xd, int yd) {
157 int x = 0, y = 1;
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);
161 }
162 }
163}
164
165
166
167
168int main(int argc, char *argv[]) {
169
170 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
171
172 // This little trick lets us print to std::cout only if
173 // a (dummy) command-line argument is provided.
174 int iprint = argc - 1;
175 Teuchos::RCP<std::ostream> outStream;
176 Teuchos::oblackholestream bhs; // outputs nothing
177 if (iprint > 0)
178 outStream = Teuchos::rcp(&std::cout, false);
179 else
180 outStream = Teuchos::rcp(&bhs, false);
181
182 // Save the format state of the original std::cout.
183 Teuchos::oblackholestream oldFormatState;
184 oldFormatState.copyfmt(std::cout);
185
186 *outStream \
187 << "===============================================================================\n" \
188 << "| |\n" \
189 << "| Unit Test (Basis_HGRAD_QUAD_Cn_FEM) |\n" \
190 << "| |\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" \
194 << "| |\n" \
195 << "| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \
196 << "| |\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" \
201 << "| |\n" \
202 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
203 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
204 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
205 << "| |\n" \
206 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
207 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
208 << "| |\n" \
209 << "===============================================================================\n"\
210 << "| TEST 1: Patch test |\n"\
211 << "===============================================================================\n";
212
213
214 int errorFlag = 0;
215
216 outStream -> precision(16);
217
218
219 try {
220
221 int max_order = 2; // max total order of polynomial solution
222 DefaultCubatureFactory<double> cubFactory; // create cubature factory
223 shards::CellTopology cell(shards::getCellTopologyData< shards::Quadrilateral<> >()); // create parent cell topology
224 shards::CellTopology side(shards::getCellTopologyData< shards::Line<> >()); // create relevant subcell (side) topology
225 int cellDim = cell.getDimension();
226 int sideDim = side.getDimension();
227
228 // Define array containing points at which the solution is evaluated, in reference cell.
229 int numIntervals = 10;
230 int numInterpPoints = (numIntervals + 1)*(numIntervals + 1);
231 FieldContainer<double> interp_points_ref(numInterpPoints, 2);
232 int counter = 0;
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;
237 counter++;
238 }
239 }
240
241 /* Parent cell definition. */
242 FieldContainer<double> cell_nodes[2];
243 cell_nodes[0].resize(1, 4, cellDim);
244 cell_nodes[1].resize(1, 4, cellDim);
245
246 // Generic parallelogram.
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;
255 // Reference quad.
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;
264
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");
268
269 for (int pcell = 0; pcell < 2; pcell++) {
270 *outStream << mystream[pcell].str();
271 FieldContainer<double> interp_points(1, numInterpPoints, cellDim);
272 CellTools<double>::mapToPhysicalFrame(interp_points, interp_points_ref, cell_nodes[pcell], cell);
273 interp_points.resize(numInterpPoints, cellDim);
274
275 for (int x_order=0; x_order <= max_order; x_order++) {
276 int max_y_order = max_order;
277 if (pcell == 0) {
278 max_y_order -= x_order;
279 }
280 for (int y_order=0; y_order <= max_y_order; y_order++) {
281
282 // evaluate exact solution
283 FieldContainer<double> exact_solution(1, numInterpPoints);
284 u_exact(exact_solution, interp_points, x_order, y_order);
285
286 int basis_order = 2;
287
288 // set test tolerance
289 double zero = basis_order*basis_order*100*INTREPID_TOL;
290
291 //create basis
292 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
293 Teuchos::rcp(new Basis_HGRAD_QUAD_Cn_FEM<double,FieldContainer<double> >(2,POINTTYPE_SPECTRAL) );
294 int numFields = basis->getCardinality();
295
296 // create cubatures
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();
301
302 /* Computational arrays. */
303 /* Section 1: Related to parent cell integration. */
304 FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim);
305 FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim);
306 FieldContainer<double> cub_weights_cell(numCubPointsCell);
307 FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim);
308 FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim);
309 FieldContainer<double> jacobian_det_cell(1, numCubPointsCell);
310 FieldContainer<double> weighted_measure_cell(1, numCubPointsCell);
311
312 FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell);
313 FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
314 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
315 FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim);
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);
318 FieldContainer<double> fe_matrix(1, numFields, numFields);
319
320 FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell);
321 FieldContainer<double> rhs_and_soln_vector(1, numFields);
322
323 /* Section 2: Related to subcell (side) integration. */
324 unsigned numSides = 4;
325 FieldContainer<double> cub_points_side(numCubPointsSide, sideDim);
326 FieldContainer<double> cub_weights_side(numCubPointsSide);
327 FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim);
328 FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim);
329 FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim);
330 FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide);
331 FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide);
332
333 FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide);
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);
336 FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide);
337 FieldContainer<double> neumann_fields_per_side(1, numFields);
338
339 /* Section 3: Related to global interpolant. */
340 FieldContainer<double> value_of_basis_at_interp_points(numFields, numInterpPoints);
341 FieldContainer<double> transformed_value_of_basis_at_interp_points(1, numFields, numInterpPoints);
342 FieldContainer<double> interpolant(1, numInterpPoints);
343
344 FieldContainer<int> ipiv(numFields);
345
346
347
348 /******************* START COMPUTATION ***********************/
349
350 // get cubature points and weights
351 cellCub->getCubature(cub_points_cell, cub_weights_cell);
352
353 // compute geometric cell information
354 CellTools<double>::setJacobian(jacobian_cell, cub_points_cell, cell_nodes[pcell], cell);
355 CellTools<double>::setJacobianInv(jacobian_inv_cell, jacobian_cell);
356 CellTools<double>::setJacobianDet(jacobian_det_cell, jacobian_cell);
357
358 // compute weighted measure
359 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
360
362 // Computing mass matrices:
363 // tabulate values of basis functions at (reference) cubature points
364 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
365
366 // transform values of basis functions
367 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
368 value_of_basis_at_cub_points_cell);
369
370 // multiply with weighted measure
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);
374
375 // compute mass matrices
377 transformed_value_of_basis_at_cub_points_cell,
378 weighted_transformed_value_of_basis_at_cub_points_cell,
379 COMP_BLAS);
381
383 // Computing stiffness matrices:
384 // tabulate gradients of basis functions at (reference) cubature points
385 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
386
387 // transform gradients of basis functions
388 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
389 jacobian_inv_cell,
390 grad_of_basis_at_cub_points_cell);
391
392 // multiply with weighted measure
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);
396
397 // compute stiffness matrices and sum into fe_matrix
399 transformed_grad_of_basis_at_cub_points_cell,
400 weighted_transformed_grad_of_basis_at_cub_points_cell,
401 COMP_BLAS,
402 true);
404
406 // Computing RHS contributions:
407 // map cell (reference) cubature points to physical space
408 CellTools<double>::mapToPhysicalFrame(cub_points_cell_physical, cub_points_cell, cell_nodes[pcell], cell);
409
410 // evaluate rhs function
411 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order);
412
413 // compute rhs
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,
417 COMP_BLAS);
418
419 // compute neumann b.c. contributions and adjust rhs
420 sideCub->getCubature(cub_points_side, cub_weights_side);
421 for (unsigned i=0; i<numSides; i++) {
422 // compute geometric cell information
423 CellTools<double>::mapToReferenceSubcell(cub_points_side_refcell, cub_points_side, sideDim, (int)i, cell);
424 CellTools<double>::setJacobian(jacobian_side_refcell, cub_points_side_refcell, cell_nodes[pcell], cell);
425 CellTools<double>::setJacobianDet(jacobian_det_side_refcell, jacobian_side_refcell);
426
427 // compute weighted edge measure
428 FunctionSpaceTools::computeEdgeMeasure<double>(weighted_measure_side_refcell,
429 jacobian_side_refcell,
430 cub_weights_side,
431 i,
432 cell);
433
434 // tabulate values of basis functions at side cubature points, in the reference parent cell domain
435 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
436 // transform
437 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
438 value_of_basis_at_cub_points_side_refcell);
439
440 // multiply with weighted measure
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);
444
445 // compute Neumann data
446 // map side cubature points in reference parent cell domain to physical space
447 CellTools<double>::mapToPhysicalFrame(cub_points_side_physical, cub_points_side_refcell, cell_nodes[pcell], cell);
448 // now compute data
449 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
450 cell, (int)i, x_order, y_order);
451
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,
455 COMP_BLAS);
456
457 // adjust RHS
458 RealSpaceTools<double>::add(rhs_and_soln_vector, neumann_fields_per_side);;
459 }
461
463 // Solution of linear system:
464 int info = 0;
465 Teuchos::LAPACK<int, double> solver;
466 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
468
470 // Building interpolant:
471 // evaluate basis at interpolation points
472 basis->getValues(value_of_basis_at_interp_points, interp_points_ref, OPERATOR_VALUE);
473 // transform values of basis functions
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);
478
479 /******************* END COMPUTATION ***********************/
480
481 RealSpaceTools<double>::subtract(interpolant, exact_solution);
482
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 << ": "
485 << RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
486 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) << "\n";
487
488 if (RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
489 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) > zero) {
490 *outStream << "\n\nPatch test failed for solution polynomial order ("
491 << x_order << ", " << y_order << ") and basis order " << basis_order << "\n\n";
492 errorFlag++;
493 }
494 } // end for y_order
495 } // end for x_order
496 } // end for pcell
497
498 }
499 // Catch unexpected errors
500 catch (const std::logic_error & err) {
501 *outStream << err.what() << "\n\n";
502 errorFlag = -1000;
503 };
504
505 if (errorFlag != 0)
506 std::cout << "End Result: TEST FAILED\n";
507 else
508 std::cout << "End Result: TEST PASSED\n";
509
510 // reset format state of std::cout
511 std::cout.copyfmt(oldFormatState);
512
513 return errorFlag;
514}
void neumann(FieldContainer< double > &, const FieldContainer< double > &, const FieldContainer< double > &, const shards::CellTopology &, int, int, int)
neumann boundary conditions
Definition test_02.cpp:112
void u_exact(FieldContainer< double > &, const FieldContainer< double > &, int, int)
exact solution
Definition test_02.cpp:156
void rhsFunc(FieldContainer< double > &, const FieldContainer< double > &, int, int)
right-hand side function
Definition test_02.cpp:76
Header file for utility class to provide array tools, such as tensor contractions,...
Header file for the Intrepid::CellTools class.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::FunctionSpaceTools class.
Header file for the Intrepid::HGRAD_QUAD_Cn_FEM class.
Header file for utility class to provide point tools, such as barycentric coordinates,...
Header file for classes providing basic linear algebra functionality in 1D, 2D and 3D.
static const double INTREPID_TOL
General purpose tolerance in, e.g., internal Newton's method to invert ref to phys maps.
static void mapToReferenceSubcell(ArraySubcellPoint &refSubcellPoints, const ArrayParamPoint &paramPoints, const int subcellDim, const int subcellOrd, const shards::CellTopology &parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
static void mapToPhysicalFrame(ArrayPhysPoint &physPoints, const ArrayRefPoint &refPoints, const ArrayCell &cellWorkset, const shards::CellTopology &cellTopo, const int &whichCell=-1)
Computes F, the reference-to-physical frame map.
static void getPhysicalSideNormals(ArraySideNormal &sideNormals, const ArrayJac &worksetJacobians, const int &worksetSideOrd, const shards::CellTopology &parentCell)
Computes non-normalized normal vectors to physical sides in a side workset .
static void setJacobianDet(ArrayJacDet &jacobianDet, const ArrayJac &jacobian)
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F.
static void setJacobianInv(ArrayJacInv &jacobianInv, const ArrayJac &jacobian)
Computes the inverse of the Jacobian matrix DF of the reference-to-physical frame map F.
A factory class that generates specific instances of cubatures.
Teuchos::RCP< Cubature< Scalar, ArrayPoint, ArrayWeight > > create(const shards::CellTopology &cellTopology, const std::vector< int > &degree)
Factory method.
Implementation of a templated lexicographical container for a multi-indexed scalar quantity....
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0.
int dimension(const int whichDim) const
Returns the specified dimension.
static void integrate(Intrepid::FieldContainer< Scalar > &outputValues, const Intrepid::FieldContainer< Scalar > &leftValues, const Intrepid::FieldContainer< Scalar > &rightValues, const ECompEngine compEngine, const bool sumInto=false)
Contracts leftValues and rightValues arrays on the point and possibly space dimensions and stores the...
static void HGRADtransformVALUE(ArrayTypeOut &outVals, const ArrayTypeIn &inVals)
Transformation of a (scalar) value field in the H-grad space, defined at points on a reference cell,...
static void scalarMultiplyDataData(ArrayOutData &outputData, ArrayInDataLeft &inputDataLeft, ArrayInDataRight &inputDataRight, const bool reciprocal=false)
Scalar multiplication of data and data; please read the description below.
static void evaluate(ArrayOutPointVals &outPointVals, const ArrayInCoeffs &inCoeffs, const ArrayInFields &inFields)
Computes point values outPointVals of a discrete function specified by the basis inFields and coeffic...
static void dotMultiplyDataData(ArrayOutData &outputData, const ArrayInDataLeft &inputDataLeft, const ArrayInDataRight &inputDataRight)
Dot product of data and data; please read the description below.
static void multiplyMeasure(ArrayTypeOut &outVals, const ArrayTypeMeasure &inMeasure, const ArrayTypeIn &inVals)
Multiplies fields inVals by weighted measures inMeasure and returns the field array outVals; this is ...
static void HGRADtransformGRAD(ArrayTypeOut &outVals, const ArrayTypeJac &jacobianInverse, const ArrayTypeIn &inVals, const char transpose='T')
Transformation of a gradient field in the H-grad space, defined at points on a reference cell,...
static void computeCellMeasure(ArrayOut &outVals, const ArrayDet &inDet, const ArrayWeights &inWeights)
Returns the weighted integration measures outVals with dimensions (C,P) used for the computation of c...
static void computeEdgeMeasure(ArrayOut &outVals, const ArrayJac &inJac, const ArrayWeights &inWeights, const int whichEdge, const shards::CellTopology &parentCell)
Returns the weighted integration measures outVals with dimensions (C,P) used for the computation of e...
static void subtract(Scalar *diffArray, const Scalar *inArray1, const Scalar *inArray2, const int size)
Subtracts contiguous data inArray2 from inArray1 of size size: diffArray = inArray1 - inArray2.
static Scalar vectorNorm(const Scalar *inVec, const size_t dim, const ENorm normType)
Computes norm (1, 2, infinity) of the vector inVec of size dim.
static void add(Scalar *sumArray, const Scalar *inArray1, const Scalar *inArray2, const int size)
Adds contiguous data inArray1 and inArray2 of size size: sumArray = inArray1 + inArray2.