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
50#include "Teuchos_oblackholestream.hpp"
51#include "Teuchos_RCP.hpp"
52#include "Teuchos_GlobalMPISession.hpp"
53
54using namespace std;
55using namespace Intrepid;
56
57#define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \
58{ \
59 ++nException; \
60 try { \
61 S ; \
62 } \
63 catch (const std::logic_error & err) { \
64 ++throwCounter; \
65 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
66 *outStream << err.what() << '\n'; \
67 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \
68 }; \
69}
70
71int main(int argc, char *argv[]) {
72
73 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
74
75 // This little trick lets us print to std::cout only if
76 // a (dummy) command-line argument is provided.
77 int iprint = argc - 1;
78 Teuchos::RCP<std::ostream> outStream;
79 Teuchos::oblackholestream bhs; // outputs nothing
80 if (iprint > 0)
81 outStream = Teuchos::rcp(&std::cout, false);
82 else
83 outStream = Teuchos::rcp(&bhs, false);
84
85 // Save the format state of the original std::cout.
86 Teuchos::oblackholestream oldFormatState;
87 oldFormatState.copyfmt(std::cout);
88
89 *outStream \
90 << "===============================================================================\n" \
91 << "| |\n" \
92 << "| Unit Test (Basis_HCURL_HEX_I1_FEM) |\n" \
93 << "| |\n" \
94 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
95 << "| 2) Basis values for VALUE and CURL operators |\n" \
96 << "| |\n" \
97 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
98 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
99 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
100 << "| |\n" \
101 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
102 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
103 << "| |\n" \
104 << "===============================================================================\n"\
105 << "| TEST 1: Basis creation, exception testing |\n"\
106 << "===============================================================================\n";
107
108 // Define basis and error flag
110 int errorFlag = 0;
111
112 // Initialize throw counter for exception testing
113 int nException = 0;
114 int throwCounter = 0;
115
116 // Define array containing the 8 vertices of the reference HEX, its center and 6 face centers
117 FieldContainer<double> hexNodes(15, 3);
118 hexNodes(0,0) = -1.0; hexNodes(0,1) = -1.0; hexNodes(0,2) = -1.0;
119 hexNodes(1,0) = 1.0; hexNodes(1,1) = -1.0; hexNodes(1,2) = -1.0;
120 hexNodes(2,0) = 1.0; hexNodes(2,1) = 1.0; hexNodes(2,2) = -1.0;
121 hexNodes(3,0) = -1.0; hexNodes(3,1) = 1.0; hexNodes(3,2) = -1.0;
122
123 hexNodes(4,0) = -1.0; hexNodes(4,1) = -1.0; hexNodes(4,2) = 1.0;
124 hexNodes(5,0) = 1.0; hexNodes(5,1) = -1.0; hexNodes(5,2) = 1.0;
125 hexNodes(6,0) = 1.0; hexNodes(6,1) = 1.0; hexNodes(6,2) = 1.0;
126 hexNodes(7,0) = -1.0; hexNodes(7,1) = 1.0; hexNodes(7,2) = 1.0;
127
128 hexNodes(8,0) = 0.0; hexNodes(8,1) = 0.0; hexNodes(8,2) = 0.0;
129
130 hexNodes(9,0) = 1.0; hexNodes(9,1) = 0.0; hexNodes(9,2) = 0.0;
131 hexNodes(10,0)= -1.0; hexNodes(10,1)= 0.0; hexNodes(10,2)= 0.0;
132
133 hexNodes(11,0)= 0.0; hexNodes(11,1)= 1.0; hexNodes(11,2)= 0.0;
134 hexNodes(12,0)= 0.0; hexNodes(12,1)= -1.0; hexNodes(12,2)= 0.0;
135
136 hexNodes(13,0)= 0.0; hexNodes(13,1)= 0.0; hexNodes(13,2)= 1.0;
137 hexNodes(14,0)= 0.0; hexNodes(14,1)= 0.0; hexNodes(14,2)= -1.0;
138
139
140 // Generic array for the output values; needs to be properly resized depending on the operator type
142
143 try{
144 // exception #1: GRAD cannot be applied to HCURL functions
145 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
146 vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0), 4 );
147 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_GRAD), throwCounter, nException );
148
149 // exception #2: DIV cannot be applied to HCURL functions
150 // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
151 vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0) );
152 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_DIV), throwCounter, nException );
153
154 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
155 // getDofTag() to access invalid array elements thereby causing bounds check exception
156 // exception #3
157 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(3,0,0), throwCounter, nException );
158 // exception #4
159 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(1,1,1), throwCounter, nException );
160 // exception #5
161 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(0,4,1), throwCounter, nException );
162 // exception #6
163 INTREPID_TEST_COMMAND( hexBasis.getDofTag(12), throwCounter, nException );
164 // exception #7
165 INTREPID_TEST_COMMAND( hexBasis.getDofTag(-1), throwCounter, nException );
166
167#ifdef HAVE_INTREPID_DEBUG
168 // Exceptions 8-15 test exception handling with incorrectly dimensioned input/output arrays
169 // exception #8: input points array must be of rank-2
170 FieldContainer<double> badPoints1(4, 5, 3);
171 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
172
173 // exception #9 dimension 1 in the input point array must equal space dimension of the cell
174 FieldContainer<double> badPoints2(4, 2);
175 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
176
177 // exception #10 output values must be of rank-3 for OPERATOR_VALUE
178 FieldContainer<double> badVals1(4, 3);
179 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_VALUE), throwCounter, nException );
180
181 // exception #11 output values must be of rank-3 for OPERATOR_CURL
182 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_CURL), throwCounter, nException );
183
184 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions)
185 FieldContainer<double> badVals2(hexBasis.getCardinality() + 1, hexNodes.dimension(0), 3);
186 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_VALUE), throwCounter, nException ) ;
187
188 // exception #13 incorrect 1st dimension of output array (must equal number of points)
189 FieldContainer<double> badVals3(hexBasis.getCardinality(), hexNodes.dimension(0) + 1, 3);
190 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals3, hexNodes, OPERATOR_VALUE), throwCounter, nException ) ;
191
192 // exception #14: incorrect 2nd dimension of output array (must equal the space dimension)
193 FieldContainer<double> badVals4(hexBasis.getCardinality(), hexNodes.dimension(0), 4);
194 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_VALUE), throwCounter, nException ) ;
195
196 // exception #15: incorrect 2nd dimension of output array (must equal the space dimension)
197 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_CURL), throwCounter, nException ) ;
198
199 // exception #16: D2 cannot be applied to HCURL functions
200 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
201 vals.resize(hexBasis.getCardinality(),
202 hexNodes.dimension(0),
203 Intrepid::getDkCardinality(OPERATOR_D2, hexBasis.getBaseCellTopology().getDimension()));
204 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_D2), throwCounter, nException );
205
206#endif
207
208 }
209 catch (const std::logic_error & err) {
210 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
211 *outStream << err.what() << '\n';
212 *outStream << "-------------------------------------------------------------------------------" << "\n\n";
213 errorFlag = -1000;
214 };
215
216 // Check if number of thrown exceptions matches the one we expect
217 // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match.
218 if (throwCounter != nException) {
219 errorFlag++;
220 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
221 }
222//#endif
223
224 *outStream \
225 << "\n"
226 << "===============================================================================\n"\
227 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
228 << "===============================================================================\n";
229
230 try{
231 std::vector<std::vector<int> > allTags = hexBasis.getAllDofTags();
232
233 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
234 for (unsigned i = 0; i < allTags.size(); i++) {
235 int bfOrd = hexBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
236
237 std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
238 if( !( (myTag[0] == allTags[i][0]) &&
239 (myTag[1] == allTags[i][1]) &&
240 (myTag[2] == allTags[i][2]) &&
241 (myTag[3] == allTags[i][3]) ) ) {
242 errorFlag++;
243 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
244 *outStream << " getDofOrdinal( {"
245 << allTags[i][0] << ", "
246 << allTags[i][1] << ", "
247 << allTags[i][2] << ", "
248 << allTags[i][3] << "}) = " << bfOrd <<" but \n";
249 *outStream << " getDofTag(" << bfOrd << ") = { "
250 << myTag[0] << ", "
251 << myTag[1] << ", "
252 << myTag[2] << ", "
253 << myTag[3] << "}\n";
254 }
255 }
256
257 // Now do the same but loop over basis functions
258 for( int bfOrd = 0; bfOrd < hexBasis.getCardinality(); bfOrd++) {
259 std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
260 int myBfOrd = hexBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
261 if( bfOrd != myBfOrd) {
262 errorFlag++;
263 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
264 *outStream << " getDofTag(" << bfOrd << ") = { "
265 << myTag[0] << ", "
266 << myTag[1] << ", "
267 << myTag[2] << ", "
268 << myTag[3] << "} but getDofOrdinal({"
269 << myTag[0] << ", "
270 << myTag[1] << ", "
271 << myTag[2] << ", "
272 << myTag[3] << "} ) = " << myBfOrd << "\n";
273 }
274 }
275 }
276 catch (const std::logic_error & err){
277 *outStream << err.what() << "\n\n";
278 errorFlag = -1000;
279 };
280
281 *outStream \
282 << "\n"
283 << "===============================================================================\n"\
284 << "| TEST 3: correctness of basis function values |\n"\
285 << "===============================================================================\n";
286
287 outStream -> precision(20);
288
289 // VALUE: Each row pair gives the 12x3 correct basis set values at an evaluation point: (P,F,D) layout
290 double basisValues[] = {
291 // bottom 4 vertices
292 0.5,0.,0., 0.,0.,0., 0.,0.,0., 0.,-0.5,0., 0.,0.,0., 0.,0.,0.,
293 0.,0.,0., 0.,0.,0., 0.,0.,0.5, 0.,0.,0., 0.,0.,0., 0.,0.,0.,
294
295 0.5,0.,0., 0.,0.5,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.,
296 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.5, 0.,0.,0., 0.,0.,0.,
297
298 0.,0.,0., 0.,0.5,0., -0.5,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.,
299 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.5, 0.,0.,0.,
300
301 0.,0.,0., 0.,0.,0., -0.5,0.,0., 0.,-0.5,0., 0.,0.,0., 0.,0.,0.,
302 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.5,
303
304 // top 4 vertices
305 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.5,0.,0., 0.,0.,0.,
306 0.,0.,0., 0.,-0.5,0., 0.,0.,0.5, 0.,0.,0., 0.,0.,0., 0.,0.,0.,
307
308 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.5,0.,0., 0.,0.5,0.,
309 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.5, 0.,0.,0., 0.,0.,0.,
310
311 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.5,0.,
312 -0.5,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.5, 0.,0.,0.,
313
314 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.,
315 -0.5,0.,0., 0.,-0.5,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.5,
316
317 // center {0, 0, 0}
318 0.125,0.,0., 0.,0.125,0., -0.125,0.,0., 0.,-0.125,0., 0.125,0.,0., 0.,0.125,0.,
319 -0.125,0.,0., 0.,-0.125,0., 0.,0.,0.125, 0.,0.,0.125, 0.,0.,0.125, 0.,0.,0.125,
320
321 // faces { 1, 0, 0} and {-1, 0, 0}
322 0.125,0.,0., 0.,0.25,0., -0.125,0.,0., 0.,0.,0., 0.125,0.,0., 0.,0.25,0.,
323 -0.125,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.25, 0.,0.,0.25, 0.,0.,0.,
324
325 0.125,0.,0., 0.,0.,0., -0.125,0.,0., 0.,-0.25,0., 0.125,0.,0., 0.,0.,0.,
326 -0.125,0.,0., 0.,-0.25,0., 0.,0.,0.25, 0.,0.,0., 0.,0.,0., 0.,0.,0.25,
327
328 // faces { 0, 1, 0} and { 0,-1, 0}
329 0.,0.,0., 0.,0.125,0., -0.25,0.,0., 0.,-0.125,0., 0.,0.,0., 0.,0.125,0.,
330 -0.25,0.,0., 0.,-0.125,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.25, 0.,0.,0.25,
331
332 0.25,0.,0., 0.,0.125,0., 0.,0.,0., 0.,-0.125,0., 0.25,0.,0., 0.,0.125,0.,
333 0.,0.,0., 0.,-0.125,0., 0.,0.,0.25, 0.,0.,0.25, 0.,0.,0., 0.,0.,0.,
334
335 // faces {0, 0, 1} and {0, 0, -1}
336 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.25,0.,0., 0.,0.25,0.,
337 -0.25,0.,0., 0.,-0.25,0., 0.,0.,0.125, 0.,0.,0.125, 0.,0.,0.125, 0.,0.,0.125,
338
339 0.25,0.,0., 0.,0.25,0., -0.25,0.,0., 0.,-0.25,0., 0.,0.,0., 0.,0.,0.,
340 0.0,0.,0., 0.,0.,0.0, 0.,0.,0.125, 0.,0.,0.125, 0.,0.,0.125, 0.,0.,0.125
341 };
342
343 // CURL: each row pair gives the 3x12 correct values of the curls of the 12 basis functions: (P,F,D) layout
344 double basisCurls[] = {
345 // bottom 4 vertices
346 0.,-0.25,0.25, 0.,0.,0.25, 0.,0.,0.25, -0.25,0.,0.25, 0.,0.25,0., 0.,0.,0.,
347 0.,0.,0., 0.25,0.,0., -0.25,0.25,0., 0.,-0.25,0., 0.,0.,0., 0.25,0.,0.,
348
349 0.,-0.25,0.25, 0.25,0.,0.25, 0.,0.,0.25, 0.,0.,0.25, 0.,0.25,0., -0.25,0.,0.,
350 0.,0.,0., 0.,0.,0., 0.,0.25,0., -0.25,-0.25,0., 0.25,0.,0., 0.,0.,0.,
351
352 0.,0.,0.25, 0.25,0.,0.25, 0.,0.25,0.25, 0.,0.,0.25, 0.,0.,0., -0.25,0.,0.,
353 0.,-0.25,0., 0.,0.,0., 0.,0.,0., -0.25,0.,0., 0.25,-0.25,0., 0.,0.25,0.,
354
355 0.,0.,0.25, 0.,0.,0.25, 0.,0.25,0.25, -0.25,0.,0.25, 0.,0.,0., 0.,0.,0.,
356 0.,-0.25,0., 0.25,0.,0., -0.25,0.,0., 0.,0.,0., 0.,-0.25,0., 0.25,0.25,0.,
357
358 // top 4 vertices
359 0.,-0.25,0., 0.,0.,0., 0.,0.,0., -0.25,0.,0., 0.,0.25,0.25, 0.,0.,0.25,
360 0.,0.,0.25, 0.25,0.,0.25, -0.25,0.25,0., 0.,-0.25,0., 0.,0.,0., 0.25,0.,0.,
361
362 0.,-0.25,0., 0.25,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.25,0.25, -0.25,0.,0.25,
363 0.,0.,0.25, 0.,0.,0.25, 0.,0.25,0., -0.25,-0.25,0., 0.25,0.,0., 0.,0.,0.,
364
365 0.,0.,0., 0.25,0.,0., 0.,0.25,0., 0.,0.,0., 0.,0.,0.25, -0.25,0.,0.25,
366 0.,-0.25,0.25, 0.,0.,0.25, 0.,0.,0., -0.25,0.,0., 0.25,-0.25,0., 0.,0.25,0.,
367
368 0.,0.,0., 0.,0.,0., 0.,0.25,0., -0.25,0.,0., 0.,0.,0.25, 0.,0.,0.25,
369 0.,-0.25,0.25, 0.25,0.,0.25, -0.25,0.,0., 0.,0.,0., 0.,-0.25,0., 0.25,0.25,0.,
370
371 // center {0, 0, 0}
372 0.,-0.125,0.125, 0.125,0.,0.125, 0.,0.125,0.125, -0.125,0.,0.125, 0.,0.125,0.125, -0.125,0.,0.125,
373 0.,-0.125,0.125, 0.125,0.,0.125, -0.125,0.125,0., -0.125,-0.125,0., 0.125,-0.125,0., 0.125,0.125,0.,
374
375 // faces { 1, 0, 0} and {-1, 0, 0}
376 0.,-0.125,0.125, 0.25,0.,0.125, 0.,0.125,0.125, 0.,0.,0.125, 0.,0.125,0.125, -0.25,0.,0.125,
377 0.,-0.125,0.125, 0.,0.,0.125, 0.,0.125,0., -0.25,-0.125,0., 0.25,-0.125,0., 0.,0.125,0.,
378
379 0.,-0.125,0.125, 0.,0.,0.125, 0.,0.125,0.125, -0.25,0.,0.125, 0.,0.125,0.125, 0.,0.,0.125,
380 0.,-0.125,0.125, 0.25,0.,0.125, -0.25,0.125,0., 0.,-0.125,0., 0.,-0.125,0., 0.25,0.125,0.,
381
382 // faces { 0, 1, 0} and { 0,-1, 0}
383 0.,0.,0.125, 0.125,0.,0.125, 0.,0.25,0.125, -0.125,0.,0.125, 0.,0.,0.125, -0.125,0.,0.125,
384 0.,-0.25,0.125, 0.125,0.,0.125, -0.125,0.,0., -0.125,0.,0., 0.125,-0.25,0., 0.125,0.25,0.,
385
386 0.,-0.25,0.125, 0.125,0.,0.125, 0.,0.,0.125, -0.125,0.,0.125, 0.,0.25,0.125, -0.125,0.,0.125,
387 0.,0.,0.125, 0.125,0.,0.125, -0.125,0.25,0., -0.125,-0.25,0., 0.125,0.,0., 0.125,0.,0.,
388
389 // faces {0, 0, 1} and {0, 0, -1}
390 0.,-0.125,0., 0.125,0.,0., 0.,0.125,0., -0.125,0.,0., 0.,0.125,0.25, -0.125,0.,0.25,
391 0.,-0.125,0.25, 0.125,0.,0.25, -0.125,0.125,0., -0.125,-0.125,0., 0.125,-0.125,0., 0.125,0.125,0.,
392
393 0.,-0.125,0.25, 0.125,0.,0.25, 0.,0.125,0.25, -0.125,0.,0.25, 0.,0.125,0., -0.125,0.,0.,
394 0.,-0.125,0., 0.125,0.,0., -0.125,0.125,0., -0.125,-0.125,0., 0.125,-0.125,0., 0.125,0.125,0.
395 };
396
397
398 try{
399
400 // Dimensions for the output arrays:
401 int numFields = hexBasis.getCardinality();
402 int numPoints = hexNodes.dimension(0);
403 int spaceDim = hexBasis.getBaseCellTopology().getDimension();
404
405 // Generic array for values and curls that will be properly sized before each call
407
408 // Check VALUE of basis functions: resize vals to rank-3 container:
409 vals.resize(numFields, numPoints, spaceDim);
410 hexBasis.getValues(vals, hexNodes, OPERATOR_VALUE);
411 for (int i = 0; i < numFields; i++) {
412 for (int j = 0; j < numPoints; j++) {
413 for (int k = 0; k < spaceDim; k++) {
414
415 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k
416 int l = k + i * spaceDim + j * spaceDim * numFields;
417 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
418 errorFlag++;
419 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
420
421 // Output the multi-index of the value where the error is:
422 *outStream << " At multi-index { ";
423 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
424 *outStream << "} computed value: " << vals(i,j,k)
425 << " but reference value: " << basisValues[l] << "\n";
426 }
427 }
428 }
429 }
430
431 // Check CURL of basis function: resize vals to rank-3 container
432 vals.resize(numFields, numPoints, spaceDim);
433 hexBasis.getValues(vals, hexNodes, OPERATOR_CURL);
434 for (int i = 0; i < numFields; i++) {
435 for (int j = 0; j < numPoints; j++) {
436 for (int k = 0; k < spaceDim; k++) {
437
438 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k
439 int l = k + i * spaceDim + j * spaceDim * numFields;
440 if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) {
441 errorFlag++;
442 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
443
444 // Output the multi-index of the value where the error is:
445 *outStream << " At multi-index { ";
446 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
447 *outStream << "} computed curl component: " << vals(i,j,k)
448 << " but reference curl component: " << basisCurls[l] << "\n";
449 }
450 }
451 }
452 }
453
454 }
455
456 // Catch unexpected errors
457 catch (const std::logic_error & err) {
458 *outStream << err.what() << "\n\n";
459 errorFlag = -1000;
460 };
461
462 *outStream \
463 << "\n"
464 << "===============================================================================\n"\
465 << "| TEST 4: correctness of DoF locations |\n"\
466 << "===============================================================================\n";
467
468 try{
469 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis =
470 Teuchos::rcp(new Basis_HCURL_HEX_I1_FEM<double, FieldContainer<double> >);
471 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface =
472 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis);
473
474 int spaceDim = 3;
476 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality(),spaceDim); // last dimension is spatial dim
477
478 // Check exceptions.
479#ifdef HAVE_INTREPID_DEBUG
480 cvals.resize(1,2,3);
481 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
482 cvals.resize(3,2);
483 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
484 cvals.resize(4,2);
485 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
486#endif
487 cvals.resize(12,spaceDim);
488 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--;
489 // Check if number of thrown exceptions matches the one we expect
490 if (throwCounter != nException) {
491 errorFlag++;
492 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
493 }
494
495 // Check mathematical correctness
496 FieldContainer<double> tangents(basis->getCardinality(),spaceDim); // tangents at each point basis point
497 tangents(0,0) = 2.0; tangents(0,1) = 0.0; tangents(0,2) = 0.0;
498 tangents(1,0) = 0.0; tangents(1,1) = 2.0; tangents(1,2) = 0.0;
499 tangents(2,0) = -2.0; tangents(2,1) = 0.0; tangents(2,2) = 0.0;
500 tangents(3,0) = 0.0; tangents(3,1) = -2.0; tangents(3,2) = 0.0;
501 tangents(4,0) = 2.0; tangents(4,1) = 0.0; tangents(4,2) = 0.0;
502 tangents(5,0) = 0.0; tangents(5,1) = 2.0; tangents(5,2) = 0.0;
503 tangents(6,0) = -2.0; tangents(6,1) = 0.0; tangents(6,2) = 0.0;
504 tangents(7,0) = 0.0; tangents(7,1) = -2.0; tangents(7,2) = 0.0;
505 tangents(8,0) = 0.0; tangents(8,1) = 0.0; tangents(8,2) = 2.0;
506 tangents(9,0) = 0.0; tangents(9,1) = 0.0; tangents(9,2) = 2.0;
507 tangents(10,0) = 0.0; tangents(10,1) = 0.0; tangents(10,2) = 2.0;
508 tangents(11,0) = 0.0; tangents(11,1) = 0.0; tangents(11,2) = 2.0;
509
510 basis->getValues(bvals, cvals, OPERATOR_VALUE);
511 char buffer[120];
512 for (int i=0; i<bvals.dimension(0); i++) {
513 for (int j=0; j<bvals.dimension(1); j++) {
514
515 double tangent = 0.0;
516 for(int d=0;d<spaceDim;d++)
517 tangent += bvals(i,j,d)*tangents(j,d);
518
519 if ((i != j) && (std::abs(tangent - 0.0) > INTREPID_TOL)) {
520 errorFlag++;
521 sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), cvals(i,2), tangent, 0.0);
522 *outStream << buffer;
523 }
524 else if ((i == j) && (std::abs(tangent - 1.0) > INTREPID_TOL)) {
525 errorFlag++;
526 sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), cvals(i,2), tangent, 1.0);
527 *outStream << buffer;
528 }
529 }
530 }
531
532 }
533 catch (const std::logic_error & err){
534 *outStream << err.what() << "\n\n";
535 errorFlag = -1000;
536 };
537
538 if (errorFlag != 0)
539 std::cout << "End Result: TEST FAILED\n";
540 else
541 std::cout << "End Result: TEST PASSED\n";
542
543 // reset format state of std::cout
544 std::cout.copyfmt(oldFormatState);
545
546 return errorFlag;
547}
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HCURL_HEX_I1_FEM class.
static const double INTREPID_TOL
General purpose tolerance in, e.g., internal Newton's method to invert ref to phys maps.
int getDkCardinality(const EOperator operatorType, const int spaceDim)
Returns cardinality of Dk, i.e., the number of all derivatives of order k.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Hexahedron cell.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Hexahedron cell.
virtual int getDofOrdinal(const int subcDim, const int subcOrd, const int subcDofOrd)
DoF tag to ordinal lookup.
virtual const std::vector< std::vector< int > > & getAllDofTags()
Retrieves all DoF tags.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
virtual const std::vector< int > & getDofTag(const int dofOrd)
DoF ordinal to DoF tag lookup.
virtual int getCardinality() const
Returns cardinality of the basis.
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.