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"
56#include "Shards_CellTopology.hpp"
57
58#include <iostream>
59using namespace Intrepid;
60
61
67int main(int argc, char *argv[]) {
68
69 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
70
71 // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
72 int iprint = argc - 1;
73
74 Teuchos::RCP<std::ostream> outStream;
75 Teuchos::oblackholestream bhs; // outputs nothing
76
77 if (iprint > 0)
78 outStream = Teuchos::rcp(&std::cout, false);
79 else
80 outStream = Teuchos::rcp(&bhs, false);
81
82 // Save the format state of the original std::cout.
83 Teuchos::oblackholestream oldFormatState;
84 oldFormatState.copyfmt(std::cout);
85
86 *outStream \
87 << "===============================================================================\n" \
88 << "| |\n" \
89 << "| Unit Test HGRAD_TET_Cn_FEM |\n" \
90 << "| |\n" \
91 << "| 1) Tests triangular Lagrange basis |\n" \
92 << "| |\n" \
93 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
94 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
95 << "| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \
96 << "| |\n" \
97 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
98 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
99 << "| |\n" \
100 << "===============================================================================\n";
101
102 int errorFlag = 0;
103
104 // Let's instantiate a basis
105 try {
106 const int deg = 3;
107 Basis_HGRAD_TET_Cn_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_WARPBLEND );
108
109 // Get a lattice
110 const int np_lattice = PointTools::getLatticeSize( myBasis.getBaseCellTopology() , deg , 0 );
111 const int nbf = myBasis.getCardinality();
112 FieldContainer<double> lattice( np_lattice , 3 );
113 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
114 myBasis.getBaseCellTopology() ,
115 deg ,
116 0 ,
117 POINTTYPE_WARPBLEND );
118 FieldContainer<double> vals( nbf , np_lattice );
119
120 myBasis.getValues( vals , lattice , OPERATOR_VALUE );
121
122 // test for Kronecker property
123 for (int i=0;i<nbf;i++) {
124 for (int j=0;j<np_lattice;j++) {
125 if ( i==j && std::abs( vals(i,j) - 1.0 ) > 100.0 * INTREPID_TOL ) {
126 errorFlag++;
127 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
128 *outStream << " Basis function " << i << " does not have unit value at its node\n";
129 }
130 if ( i!=j && std::abs( vals(i,j) ) > 100.0 * INTREPID_TOL ) {
131 errorFlag++;
132 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
133 *outStream << " Basis function " << i << " does not vanish at node " << j << "\n";
134 *outStream << " Basis function value is " << vals(i,j) << "\n";
135 }
136 }
137 }
138 }
139 catch (const std::exception & err) {
140 *outStream << err.what() << "\n\n";
141 errorFlag = -1000;
142 }
143
144 try {
145 const int deg = 4;
146 Basis_HGRAD_TET_Cn_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_WARPBLEND );
147 const std::vector<std::vector<std::vector<int> > >& dofdata = myBasis.getDofOrdinalData();
148 for (unsigned d=0;d<dofdata.size();d++) {
149 std::cout << "Dimension " << d << "\n";
150 for (unsigned f=0;f<dofdata[d].size();f++) {
151 std::cout << "\tFacet number " << f << "\n";
152 std::cout << "\t\tDOFS:\n";
153 for (unsigned n=0;n<dofdata[d][f].size();n++) {
154 if ( dofdata[d][f][n] >= 0 ) {
155 std::cout << "\t\t\t" << dofdata[d][f][n] << "\n";
156 }
157 }
158 }
159 }
160 }
161 catch (const std::exception & err) {
162 *outStream << err.what() << "\n\n";
163 errorFlag = -1000;
164 }
165
166 try {
167 const int deg = 3;
168 Basis_HGRAD_TET_Cn_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_EQUISPACED );
169
170 // Get a lattice
171 const int np_lattice = PointTools::getLatticeSize( myBasis.getBaseCellTopology() , deg , 0 );
172 const int nbf = myBasis.getCardinality();
173 FieldContainer<double> lattice( np_lattice , 3 );
174 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
175 myBasis.getBaseCellTopology() ,
176 deg ,
177 0 ,
178 POINTTYPE_EQUISPACED );
179 FieldContainer<double> vals( nbf , np_lattice , 3 );
180
181 myBasis.getValues( vals , lattice , OPERATOR_GRAD );
182
183 }
184 catch (const std::exception & err) {
185 *outStream << err.what() << "\n\n";
186 errorFlag = -1000;
187 }
188
189 if (errorFlag != 0)
190 std::cout << "End Result: TEST FAILED\n";
191 else
192 std::cout << "End Result: TEST PASSED\n";
193
194 // reset format state of std::cout
195 std::cout.copyfmt(oldFormatState);
196
197 return errorFlag;
198}
int main(int argc, char *argv[])
Tests for Lagrange basis on tets. Tests Kronecker property of basis and basic execution of differenti...
Definition: test_01.cpp:67
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HGRAD_TET_Cn_FEM class.
Header file for utility class to provide point tools, such as barycentric coordinates,...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
virtual const std::vector< std::vector< std::vector< int > > > & getDofOrdinalData()
DoF tag to ordinal data structure.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
virtual int getCardinality() const
Returns cardinality of the basis.
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...