Intrepid
test_01.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
44
51#include "Teuchos_oblackholestream.hpp"
52#include "Teuchos_RCP.hpp"
53#include "Teuchos_GlobalMPISession.hpp"
54#include "Intrepid_HGRAD_TRI_Cn_FEM_ORTH.hpp"
57#include "Shards_CellTopology.hpp"
58#include <iostream>
59using namespace Intrepid;
60
66int main(int argc, char *argv[]) {
67
68 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
69
70 // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
71 int iprint = argc - 1;
72
73 Teuchos::RCP<std::ostream> outStream;
74 Teuchos::oblackholestream bhs; // outputs nothing
75
76 if (iprint > 0)
77 outStream = Teuchos::rcp(&std::cout, false);
78 else
79 outStream = Teuchos::rcp(&bhs, false);
80
81 // Save the format state of the original std::cout.
82 Teuchos::oblackholestream oldFormatState;
83 oldFormatState.copyfmt(std::cout);
84
85 *outStream \
86 << "===============================================================================\n" \
87 << "| |\n" \
88 << "| Unit Test OrthogonalBases |\n" \
89 << "| |\n" \
90 << "| 1) Tests orthogonality of triangular orthogonal basis (Dubiner) |\n" \
91 << "| |\n" \
92 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
93 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
94 << "| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \
95 << "| |\n" \
96 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
97 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
98 << "| |\n" \
99 << "===============================================================================\n";
100
101 int errorFlag = 0;
102
103 // First, get a reference quadrature rule
104 try {
106 FieldContainer<double> cubPts( myCub.getNumPoints() , 2 );
107 FieldContainer<double> cubWts( myCub.getNumPoints() );
108
109 myCub.getCubature( cubPts , cubWts );
110
111 // Tabulate the basis functions at the cubature points
112 const int deg =3;
114 const int polydim = myBasis.getCardinality();
115
116 FieldContainer<double> basisAtCubPts( polydim , myCub.getNumPoints() );
117
118 myBasis.getValues( basisAtCubPts , cubPts , OPERATOR_VALUE );
119
120 // Now let's compute the mass matrix
121 for (int i=0;i<polydim;i++) {
122 for (int j=i;j<polydim;j++) {
123 double cur = 0;
124 for (int k=0;k<myCub.getNumPoints();k++) {
125 cur += cubWts(k) * basisAtCubPts( i , k ) * basisAtCubPts( j , k );
126 }
127 if (i != j && fabs( cur ) > INTREPID_TOL) {
128 errorFlag++;
129 }
130 else if (i == j && fabs( cur ) < INTREPID_TOL ) {
131 errorFlag++;
132 }
133
134 }
135 }
136 }
137 catch (const std::exception & err) {
138 *outStream << err.what() << "\n\n";
139 errorFlag = -1000;
140 }
141
142 // compare the points against FIAT-tabulated values on a lattice
143 try {
144 const int deg = 3;
146 shards::CellTopology myTri_3( shards::getCellTopologyData< shards::Triangle<3> >() );
147 const int np_lattice = PointTools::getLatticeSize( myTri_3 , deg , 0 );
148 FieldContainer<double> lattice( np_lattice , 2);
149 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
150 myTri_3 ,
151 deg ,
152 0 ,
153 POINTTYPE_EQUISPACED );
154 const int polydim = myBasis.getCardinality();
155
156 FieldContainer<double> dBasisAtLattice( polydim , np_lattice , 2 );
157 myBasis.getValues( dBasisAtLattice , lattice , OPERATOR_D1 );
158
159 double fiat_vals[] = {
1600.000000000000000e+00, 0.000000000000000e+00,
1610.000000000000000e+00, 0.000000000000000e+00,
1620.000000000000000e+00, 0.000000000000000e+00,
1630.000000000000000e+00, 0.000000000000000e+00,
1640.000000000000000e+00, 0.000000000000000e+00,
1650.000000000000000e+00, 0.000000000000000e+00,
1660.000000000000000e+00, 0.000000000000000e+00,
1670.000000000000000e+00, 0.000000000000000e+00,
1680.000000000000000e+00, 0.000000000000000e+00,
1690.000000000000000e+00, 0.000000000000000e+00,
1703.464101615137754e+00, 1.732050807568877e+00,
1713.464101615137754e+00, 1.732050807568877e+00,
1723.464101615137754e+00, 1.732050807568877e+00,
1733.464101615137754e+00, 1.732050807568877e+00,
1743.464101615137754e+00, 1.732050807568877e+00,
1753.464101615137754e+00, 1.732050807568877e+00,
1763.464101615137754e+00, 1.732050807568877e+00,
1773.464101615137754e+00, 1.732050807568877e+00,
1783.464101615137754e+00, 1.732050807568877e+00,
1793.464101615137754e+00, 1.732050807568877e+00,
1800.000000000000000e+00, 3.000000000000000e+00,
1810.000000000000000e+00, 3.000000000000000e+00,
1820.000000000000000e+00, 3.000000000000000e+00,
1830.000000000000000e+00, 3.000000000000000e+00,
1840.000000000000000e+00, 3.000000000000000e+00,
1850.000000000000000e+00, 3.000000000000000e+00,
1860.000000000000000e+00, 3.000000000000000e+00,
1870.000000000000000e+00, 3.000000000000000e+00,
1880.000000000000000e+00, 3.000000000000000e+00,
1890.000000000000000e+00, 3.000000000000000e+00,
190-1.643167672515498e+01, -5.477225575051661e+00,
191-5.477225575051661e+00, 0.000000000000000e+00,
1925.477225575051660e+00, 5.477225575051660e+00,
1931.643167672515498e+01, 1.095445115010332e+01,
194-1.095445115010332e+01, -3.651483716701107e+00,
195-9.121412916732176e-16, 1.825741858350553e+00,
1961.095445115010332e+01, 7.302967433402213e+00,
197-5.477225575051661e+00, -1.825741858350554e+00,
1985.477225575051660e+00, 3.651483716701107e+00,
1990.000000000000000e+00, 0.000000000000000e+00,
200-4.242640687119285e+00, -1.272792206135786e+01,
201-4.242640687119285e+00, -5.656854249492381e+00,
202-4.242640687119285e+00, 1.414213562373094e+00,
203-4.242640687119285e+00, 8.485281374238570e+00,
2042.828427124746189e+00, -5.656854249492381e+00,
2052.828427124746189e+00, 1.414213562373094e+00,
2062.828427124746189e+00, 8.485281374238568e+00,
2079.899494936611664e+00, 1.414213562373094e+00,
2089.899494936611664e+00, 8.485281374238568e+00,
2091.697056274847714e+01, 8.485281374238570e+00,
2100.000000000000000e+00, -9.797958971132712e+00,
2110.000000000000000e+00, -9.797958971132712e+00,
2120.000000000000000e+00, -9.797958971132710e+00,
2130.000000000000000e+00, -9.797958971132712e+00,
2140.000000000000000e+00, -1.632993161855452e+00,
2150.000000000000000e+00, -1.632993161855452e+00,
2160.000000000000000e+00, -1.632993161855452e+00,
2170.000000000000000e+00, 6.531972647421806e+00,
2180.000000000000000e+00, 6.531972647421806e+00,
2190.000000000000000e+00, 1.469693845669907e+01,
2204.489988864128730e+01, 1.122497216032182e+01,
221-4.988876515698587e+00, -6.236095644623235e+00,
222-4.988876515698591e+00, 1.247219128924645e+00,
2234.489988864128730e+01, 3.367491648096547e+01,
2241.995550606279436e+01, 4.988876515698590e+00,
225-4.988876515698589e+00, -2.494438257849295e+00,
2261.995550606279435e+01, 1.496662954709576e+01,
2274.988876515698590e+00, 1.247219128924648e+00,
2284.988876515698586e+00, 3.741657386773940e+00,
2290.000000000000000e+00, 0.000000000000000e+00,
2301.897366596101028e+01, 2.846049894151541e+01,
2316.324555320336759e+00, -7.378647873726218e+00,
232-6.324555320336757e+00, -1.370320319406298e+01,
233-1.897366596101028e+01, 9.486832980505138e+00,
234-1.686548085423136e+01, 4.216370213557840e+00,
235-1.404333387430680e-15, -2.108185106778921e+00,
2361.686548085423135e+01, 2.108185106778919e+01,
237-2.319003617456811e+01, -5.270462766947300e+00,
2382.319003617456811e+01, 1.791957340762081e+01,
2390.000000000000000e+00, 0.000000000000000e+00,
2404.898979485566356e+00, 3.184336665618131e+01,
2414.898979485566356e+00, 1.224744871391589e+01,
2424.898979485566356e+00, -7.348469228349532e+00,
2434.898979485566356e+00, -2.694438717061496e+01,
244-3.265986323710904e+00, -4.898979485566356e+00,
245-3.265986323710904e+00, -1.632993161855452e+00,
246-3.265986323710904e+00, 1.632993161855451e+00,
2471.143095213298816e+01, -7.348469228349537e+00,
2481.143095213298816e+01, 1.877942136133769e+01,
2494.898979485566356e+01, 2.449489742783178e+01,
2500.000000000000000e+00, 2.121320343559643e+01,
2510.000000000000000e+00, 2.121320343559643e+01,
2520.000000000000000e+00, 2.121320343559642e+01,
2530.000000000000000e+00, 2.121320343559643e+01,
2540.000000000000000e+00, -4.714045207910317e+00,
2550.000000000000000e+00, -4.714045207910317e+00,
2560.000000000000000e+00, -4.714045207910317e+00,
2570.000000000000000e+00, 2.357022603955157e+00,
2580.000000000000000e+00, 2.357022603955157e+00,
2590.000000000000000e+00, 4.242640687119285e+01
260 };
261
262 int fiat_index_cur = 0;
263 for (int i=0;i<polydim;i++) {
264 for (int j=0;j<np_lattice;j++) {
265 for (int k=0;k<2;k++) {
266 if (std::abs( dBasisAtLattice(i,j,k) - fiat_vals[fiat_index_cur] ) > INTREPID_TOL ) {
267 errorFlag++;
268 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
269
270 // Output the multi-index of the value where the error is:
271 *outStream << " At multi-index { ";
272 *outStream << i << " " << j << " " << k;
273 *outStream << "} computed value: " << dBasisAtLattice(i,j,k)
274 << " but correct value: " << fiat_vals[fiat_index_cur] << "\n";
275 *outStream << "Difference: " << std::abs( dBasisAtLattice(i,j,k) - fiat_vals[fiat_index_cur] ) << "\n";
276 }
277 fiat_index_cur++;
278 }
279 }
280 }
281 }
282 catch (const std::exception & err) {
283 *outStream << err.what() << "\n\n";
284 errorFlag = -1000;
285 }
286
287 // do second order derivatives
288 try {
289 const int deg = 3;
290 shards::CellTopology myTri_3( shards::getCellTopologyData< shards::Triangle<3> >() );
292 const int polydim = myBasis.getCardinality();
293 const int np_lattice = PointTools::getLatticeSize( myTri_3 , deg , 0 );
294 FieldContainer<double> lattice( np_lattice , 2);
295 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
296 myTri_3 ,
297 deg ,
298 0 ,
299 POINTTYPE_EQUISPACED );
300
301
302 FieldContainer<double> dBasisAtLattice( polydim , np_lattice , 3 );
303 myBasis.getValues( dBasisAtLattice , lattice , OPERATOR_D2 );
304
305 const double fiat_vals[] = {
3060.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3070.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3080.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3090.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3100.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3110.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3120.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3130.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3140.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3150.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3160.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3170.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3180.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3190.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3200.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3210.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3220.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3230.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3240.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3250.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3260.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3270.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3280.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3290.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3300.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3310.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3320.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3330.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3340.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3350.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
3363.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3373.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3383.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3393.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3403.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3413.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3423.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3433.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3443.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3453.286335345030997e+01, 1.643167672515498e+01, 5.477225575051661e+00,
3460.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3470.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3480.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3490.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3500.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3510.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3520.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3530.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3540.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3550.000000000000000e+00, 2.121320343559642e+01, 2.121320343559642e+01,
3560.000000000000000e+00, 0.000000000000000e+00, 2.449489742783178e+01,
3570.000000000000000e+00, 0.000000000000000e+00, 2.449489742783178e+01,
3580.000000000000000e+00, 0.000000000000000e+00, 2.449489742783177e+01,
3590.000000000000000e+00, 0.000000000000000e+00, 2.449489742783178e+01,
3600.000000000000000e+00, 0.000000000000000e+00, 2.449489742783178e+01,
3610.000000000000000e+00, 0.000000000000000e+00, 2.449489742783178e+01,
3620.000000000000000e+00, 0.000000000000000e+00, 2.449489742783178e+01,
3630.000000000000000e+00, 0.000000000000000e+00, 2.449489742783177e+01,
3640.000000000000000e+00, 0.000000000000000e+00, 2.449489742783177e+01,
3650.000000000000000e+00, 0.000000000000000e+00, 2.449489742783178e+01,
366-2.244994432064365e+02, -8.979977728257460e+01, -2.244994432064365e+01,
367-7.483314773547885e+01, -1.496662954709577e+01, 7.483314773547880e+00,
3687.483314773547882e+01, 5.986651818838305e+01, 3.741657386773942e+01,
3692.244994432064365e+02, 1.346996659238619e+02, 6.734983296193094e+01,
370-1.496662954709577e+02, -5.986651818838306e+01, -1.496662954709577e+01,
371-1.246222254316567e-14, 1.496662954709576e+01, 1.496662954709576e+01,
3721.496662954709576e+02, 8.979977728257458e+01, 4.489988864128730e+01,
373-7.483314773547885e+01, -2.993325909419154e+01, -7.483314773547884e+00,
3747.483314773547882e+01, 4.489988864128729e+01, 2.244994432064365e+01,
3750.000000000000000e+00, 0.000000000000000e+00, 0.000000000000000e+00,
376-3.794733192202055e+01, -1.517893276880822e+02, -9.486832980505139e+01,
377-3.794733192202055e+01, -6.324555320336759e+01, -6.324555320336759e+00,
378-3.794733192202055e+01, 2.529822128134703e+01, 8.221921916437785e+01,
379-3.794733192202055e+01, 1.138419957660617e+02, 1.707629936490925e+02,
3805.059644256269407e+01, -6.324555320336759e+01, -5.059644256269407e+01,
3815.059644256269407e+01, 2.529822128134703e+01, 3.794733192202055e+01,
3825.059644256269407e+01, 1.138419957660616e+02, 1.264911064067352e+02,
3831.391402170474087e+02, 2.529822128134704e+01, -6.324555320336762e+00,
3841.391402170474087e+02, 1.138419957660617e+02, 8.221921916437786e+01,
3852.276839915321233e+02, 1.138419957660617e+02, 3.794733192202055e+01,
3860.000000000000000e+00, -5.878775382679627e+01, -1.616663230236897e+02,
3870.000000000000000e+00, -5.878775382679627e+01, -9.308061022576078e+01,
3880.000000000000000e+00, -5.878775382679627e+01, -2.449489742783179e+01,
3890.000000000000000e+00, -5.878775382679627e+01, 4.409081537009720e+01,
3900.000000000000000e+00, 9.797958971132708e+00, -5.878775382679628e+01,
3910.000000000000000e+00, 9.797958971132708e+00, 9.797958971132703e+00,
3920.000000000000000e+00, 9.797958971132708e+00, 7.838367176906168e+01,
3930.000000000000000e+00, 7.838367176906168e+01, 4.409081537009718e+01,
3940.000000000000000e+00, 7.838367176906168e+01, 1.126765281680262e+02,
3950.000000000000000e+00, 1.469693845669907e+02, 1.469693845669907e+02,
3960.000000000000000e+00, 0.000000000000000e+00, -1.272792206135786e+02,
3970.000000000000000e+00, 0.000000000000000e+00, -1.272792206135786e+02,
3980.000000000000000e+00, 0.000000000000000e+00, -1.272792206135785e+02,
3990.000000000000000e+00, 0.000000000000000e+00, -1.272792206135786e+02,
4000.000000000000000e+00, 0.000000000000000e+00, -2.828427124746190e+01,
4010.000000000000000e+00, 0.000000000000000e+00, -2.828427124746190e+01,
4020.000000000000000e+00, 0.000000000000000e+00, -2.828427124746190e+01,
4030.000000000000000e+00, 0.000000000000000e+00, 7.071067811865474e+01,
4040.000000000000000e+00, 0.000000000000000e+00, 7.071067811865474e+01,
4050.000000000000000e+00, 0.000000000000000e+00, 1.697056274847714e+02
406
407 };
408 int fiat_index_cur = 0;
409 for (int i=0;i<polydim;i++) {
410 for (int j=0;j<np_lattice;j++) {
411 for (int k=0;k<3;k++) {
412 if (std::abs( dBasisAtLattice(i,j,k) - fiat_vals[fiat_index_cur] ) > 10.0*INTREPID_TOL ) {
413 errorFlag++;
414 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
415
416 // Output the multi-index of the value where the error is:
417 *outStream << " At multi-index { ";
418 *outStream << i << " " << j << " " << k;
419 *outStream << "} computed value: " << dBasisAtLattice(i,j,k)
420 << " but correct value: " << fiat_vals[fiat_index_cur] << "\n";
421 *outStream << "Difference: " << std::abs( dBasisAtLattice(i,j,k) - fiat_vals[fiat_index_cur] ) << "\n";
422 }
423 fiat_index_cur++;
424 }
425 }
426 }
427 }
428 catch (const std::exception & err) {
429 *outStream << err.what() << "\n\n";
430 errorFlag = -1000;
431 }
432
433
434 if (errorFlag != 0)
435 std::cout << "End Result: TEST FAILED\n";
436 else
437 std::cout << "End Result: TEST PASSED\n";
438
439 // reset format state of std::cout
440 std::cout.copyfmt(oldFormatState);
441
442 return errorFlag;
443}
int main(int argc, char *argv[])
Tests for orthogonal basis on tets. Tests diagonality of mass matrices and does a code comparison to ...
Definition: test_01.cpp:66
Header file for the Intrepid::CubatureDirectTriDefault class.
Header file for utility class to provide multidimensional containers.
Header file for utility class to provide point tools, such as barycentric coordinates,...
Implementation of the default H(grad)-compatible orthogonal basis (Dubiner) of arbitrary degree on tr...
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
virtual int getCardinality() const
Returns cardinality of the basis.
Defines direct integration rules on a triangle.
virtual void getCubature(ArrayPoint &cubPoints, ArrayWeight &cubWeights) const
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
virtual int getNumPoints() const
Returns the number of cubature points.
Implementation of a templated lexicographical container for a multi-indexed scalar quantity....
static int getLatticeSize(const shards::CellTopology &cellType, const int order, const int offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...