|
Intrepid
|
00001 // @HEADER 00002 // ************************************************************************ 00003 // 00004 // Intrepid Package 00005 // Copyright (2007) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Pavel Bochev (pbboche@sandia.gov) 00038 // Denis Ridzal (dridzal@sandia.gov), or 00039 // Kara Peterson (kjpeter@sandia.gov) 00040 // 00041 // ************************************************************************ 00042 // @HEADER 00043 00048 #include "Intrepid_FieldContainer.hpp" 00049 #include "Intrepid_HGRAD_HEX_I2_FEM.hpp" 00050 #include "Teuchos_oblackholestream.hpp" 00051 #include "Teuchos_RCP.hpp" 00052 #include "Teuchos_GlobalMPISession.hpp" 00053 00054 using namespace std; 00055 using namespace Intrepid; 00056 00057 #define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \ 00058 { \ 00059 ++nException; \ 00060 try { \ 00061 S ; \ 00062 } \ 00063 catch (std::logic_error err) { \ 00064 ++throwCounter; \ 00065 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \ 00066 *outStream << err.what() << '\n'; \ 00067 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \ 00068 }; \ 00069 } 00070 00071 int main(int argc, char *argv[]) { 00072 00073 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00074 00075 // This little trick lets us print to std::cout only if 00076 // a (dummy) command-line argument is provided. 00077 int iprint = argc - 1; 00078 Teuchos::RCP<std::ostream> outStream; 00079 Teuchos::oblackholestream bhs; // outputs nothing 00080 if (iprint > 0) 00081 outStream = Teuchos::rcp(&std::cout, false); 00082 else 00083 outStream = Teuchos::rcp(&bhs, false); 00084 00085 // Save the format state of the original std::cout. 00086 Teuchos::oblackholestream oldFormatState; 00087 oldFormatState.copyfmt(std::cout); 00088 00089 *outStream \ 00090 << "===============================================================================\n" \ 00091 << "| |\n" \ 00092 << "| Unit Test (Basis_HGRAD_HEX_I2_FEM) |\n" \ 00093 << "| |\n" \ 00094 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \ 00095 << "| 2) Basis values for VALUE, GRAD, and Dk operators |\n" \ 00096 << "| |\n" \ 00097 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00098 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00099 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00100 << "| |\n" \ 00101 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00102 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00103 << "| |\n" \ 00104 << "===============================================================================\n"\ 00105 << "| TEST 1: Basis creation, exception testing |\n"\ 00106 << "===============================================================================\n"; 00107 00108 // Define basis and error flag 00109 Basis_HGRAD_HEX_I2_FEM<double, FieldContainer<double> > hexBasis; 00110 int errorFlag = 0; 00111 00112 // Initialize throw counter for exception testing 00113 int nException = 0; 00114 int throwCounter = 0; 00115 00116 // Define arrayS containing the 20 nodes of hexahedron<20> topology 00117 FieldContainer<double> hexNodes(20, 3); 00118 // vertices 00119 hexNodes(0, 0) = -1.0; hexNodes(0, 1) = -1.0; hexNodes(0, 2) = -1.0; 00120 hexNodes(1, 0) = 1.0; hexNodes(1, 1) = -1.0; hexNodes(1, 2) = -1.0; 00121 hexNodes(2, 0) = 1.0; hexNodes(2, 1) = 1.0; hexNodes(2, 2) = -1.0; 00122 hexNodes(3, 0) = -1.0; hexNodes(3, 1) = 1.0; hexNodes(3, 2) = -1.0; 00123 00124 hexNodes(4, 0) = -1.0; hexNodes(4, 1) = -1.0; hexNodes(4, 2) = 1.0; 00125 hexNodes(5, 0) = 1.0; hexNodes(5, 1) = -1.0; hexNodes(5, 2) = 1.0; 00126 hexNodes(6, 0) = 1.0; hexNodes(6, 1) = 1.0; hexNodes(6, 2) = 1.0; 00127 hexNodes(7, 0) = -1.0; hexNodes(7, 1) = 1.0; hexNodes(7, 2) = 1.0; 00128 00129 // nodes on edges 00130 hexNodes(8, 0) = 0.0; hexNodes(8, 1) = -1.0; hexNodes(8, 2) = -1.0; 00131 hexNodes(9, 0) = 1.0; hexNodes(9, 1) = 0.0; hexNodes(9, 2) = -1.0; 00132 hexNodes(10,0) = 0.0; hexNodes(10,1) = 1.0; hexNodes(10,2) = -1.0; 00133 hexNodes(11,0) = -1.0; hexNodes(11,1) = 0.0; hexNodes(11,2) = -1.0; 00134 hexNodes(12,0) = -1.0; hexNodes(12,1) = -1.0; hexNodes(12,2) = 0.0; 00135 hexNodes(13,0) = 1.0; hexNodes(13,1) = -1.0; hexNodes(13,2) = 0.0; 00136 hexNodes(14,0) = 1.0; hexNodes(14,1) = 1.0; hexNodes(14,2) = 0.0; 00137 hexNodes(15,0) = -1.0; hexNodes(15,1) = 1.0; hexNodes(15,2) = 0.0; 00138 hexNodes(16,0) = 0.0; hexNodes(16,1) = -1.0; hexNodes(16,2) = 1.0; 00139 hexNodes(17,0) = 1.0; hexNodes(17,1) = 0.0; hexNodes(17,2) = 1.0; 00140 hexNodes(18,0) = 0.0; hexNodes(18,1) = 1.0; hexNodes(18,2) = 1.0; 00141 hexNodes(19,0) = -1.0; hexNodes(19,1) = 0.0; hexNodes(19,2) = 1.0; 00142 00143 // Generic array for the output values; needs to be properly resized depending on the operator type 00144 FieldContainer<double> vals; 00145 00146 try{ 00147 // exception #1: CURL cannot be applied to scalar functions in 3D 00148 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00149 vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0), 4 ); 00150 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_CURL), throwCounter, nException ); 00151 00152 // exception #2: DIV cannot be applied to scalar functions in 3D 00153 // resize vals to rank-2 container with dimensions (num. basis functions, num. points) 00154 vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0) ); 00155 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_DIV), throwCounter, nException ); 00156 00157 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00158 // getDofTag() to access invalid array elements thereby causing bounds check exception 00159 // exception #3 00160 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(3,10,0), throwCounter, nException ); 00161 // exception #4 00162 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(1,2,1), throwCounter, nException ); 00163 // exception #5 00164 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(0,4,1), throwCounter, nException ); 00165 // exception #6 00166 INTREPID_TEST_COMMAND( hexBasis.getDofTag(21), throwCounter, nException ); 00167 // exception #7 00168 INTREPID_TEST_COMMAND( hexBasis.getDofTag(-1), throwCounter, nException ); 00169 00170 #ifdef HAVE_INTREPID_DEBUG 00171 // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays 00172 // exception #8: input points array must be of rank-2 00173 FieldContainer<double> badPoints1(4, 5, 3); 00174 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00175 00176 // exception #9 dimension 1 in the input point array must equal space dimension of the cell 00177 FieldContainer<double> badPoints2(4, hexBasis.getBaseCellTopology().getDimension() - 1); 00178 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00179 00180 // exception #10 output values must be of rank-2 for OPERATOR_VALUE 00181 FieldContainer<double> badVals1(4, 3, 1); 00182 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_VALUE), throwCounter, nException ); 00183 00184 // exception #11 output values must be of rank-3 for OPERATOR_GRAD 00185 FieldContainer<double> badVals2(4, 3); 00186 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_GRAD), throwCounter, nException ); 00187 00188 // exception #12 output values must be of rank-3 for OPERATOR_D1 00189 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D1), throwCounter, nException ); 00190 00191 // exception #13 output values must be of rank-3 for OPERATOR_D2 00192 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D2), throwCounter, nException ); 00193 00194 // exception #14 incorrect 0th dimension of output array (must equal number of basis functions) 00195 FieldContainer<double> badVals3(hexBasis.getCardinality() + 1, hexNodes.dimension(0)); 00196 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals3, hexNodes, OPERATOR_VALUE), throwCounter, nException ); 00197 00198 // exception #15 incorrect 1st dimension of output array (must equal number of points) 00199 FieldContainer<double> badVals4(hexBasis.getCardinality(), hexNodes.dimension(0) + 1); 00200 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_VALUE), throwCounter, nException ); 00201 00202 // exception #16: incorrect 2nd dimension of output array (must equal the space dimension) 00203 FieldContainer<double> badVals5(hexBasis.getCardinality(), hexNodes.dimension(0), hexBasis.getBaseCellTopology().getDimension() - 1); 00204 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals5, hexNodes, OPERATOR_GRAD), throwCounter, nException ); 00205 00206 // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 3D) 00207 FieldContainer<double> badVals6(hexBasis.getCardinality(), hexNodes.dimension(0), 40); 00208 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals6, hexNodes, OPERATOR_D2), throwCounter, nException ); 00209 00210 // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 3D) 00211 FieldContainer<double> badVals7(hexBasis.getCardinality(), hexNodes.dimension(0), 50); 00212 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals7, hexNodes, OPERATOR_D3), throwCounter, nException ); 00213 #endif 00214 00215 } 00216 catch (std::logic_error err) { 00217 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00218 *outStream << err.what() << '\n'; 00219 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00220 errorFlag = -1000; 00221 }; 00222 00223 // Check if number of thrown exceptions matches the one we expect 00224 // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match. 00225 if (throwCounter != nException) { 00226 errorFlag++; 00227 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00228 } 00229 00230 *outStream \ 00231 << "\n" 00232 << "===============================================================================\n"\ 00233 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00234 << "===============================================================================\n"; 00235 00236 try{ 00237 std::vector<std::vector<int> > allTags = hexBasis.getAllDofTags(); 00238 00239 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00240 for (unsigned i = 0; i < allTags.size(); i++) { 00241 int bfOrd = hexBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00242 00243 std::vector<int> myTag = hexBasis.getDofTag(bfOrd); 00244 if( !( (myTag[0] == allTags[i][0]) && 00245 (myTag[1] == allTags[i][1]) && 00246 (myTag[2] == allTags[i][2]) && 00247 (myTag[3] == allTags[i][3]) ) ) { 00248 errorFlag++; 00249 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00250 *outStream << " getDofOrdinal( {" 00251 << allTags[i][0] << ", " 00252 << allTags[i][1] << ", " 00253 << allTags[i][2] << ", " 00254 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00255 *outStream << " getDofTag(" << bfOrd << ") = { " 00256 << myTag[0] << ", " 00257 << myTag[1] << ", " 00258 << myTag[2] << ", " 00259 << myTag[3] << "}\n"; 00260 } 00261 } 00262 00263 // Now do the same but loop over basis functions 00264 for( int bfOrd = 0; bfOrd < hexBasis.getCardinality(); bfOrd++) { 00265 std::vector<int> myTag = hexBasis.getDofTag(bfOrd); 00266 int myBfOrd = hexBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00267 if( bfOrd != myBfOrd) { 00268 errorFlag++; 00269 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00270 *outStream << " getDofTag(" << bfOrd << ") = { " 00271 << myTag[0] << ", " 00272 << myTag[1] << ", " 00273 << myTag[2] << ", " 00274 << myTag[3] << "} but getDofOrdinal({" 00275 << myTag[0] << ", " 00276 << myTag[1] << ", " 00277 << myTag[2] << ", " 00278 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00279 } 00280 } 00281 } 00282 catch (std::logic_error err){ 00283 *outStream << err.what() << "\n\n"; 00284 errorFlag = -1000; 00285 }; 00286 00287 00288 *outStream \ 00289 << "\n" 00290 << "===============================================================================\n"\ 00291 << "| TEST 3: correctness of basis function values |\n"\ 00292 << "===============================================================================\n"; 00293 00294 outStream -> precision(20); 00295 00296 // VALUE: Each row gives the 8 correct basis set values at an evaluation point 00297 double basisValues[] = { 00298 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00299 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00300 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00301 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00302 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00303 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00304 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00305 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00306 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00307 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00308 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00309 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00310 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, \ 00311 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, \ 00312 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, \ 00313 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, \ 00314 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, \ 00315 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, \ 00316 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, \ 00317 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0 }; 00318 00319 00320 // GRAD, D1, D2, D3 test values are stored in files due to their large size 00321 std::string fileName; 00322 std::ifstream dataFile; 00323 00324 00325 // GRAD and D1 values are stored in (F,P,D) format in a data file. Read file and do the test 00326 std::vector<double> basisGrads; // Flat array for the gradient values. 00327 00328 fileName = "./testdata/HEX_I2_GradVals.dat"; 00329 dataFile.open(fileName.c_str()); 00330 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error, 00331 ">>> ERROR (HGRAD_HEX_I2/test01): could not open GRAD values data file, test aborted."); 00332 while (!dataFile.eof() ){ 00333 double temp; 00334 string line; // string for one line of input file 00335 std::getline(dataFile, line); // get next line from file 00336 stringstream data_line(line); // convert to stringstream 00337 while(data_line >> temp){ // extract value from line 00338 basisGrads.push_back(temp); // push into vector 00339 } 00340 } 00341 // It turns out that just closing and then opening the ifstream variable does not reset it 00342 // and subsequent open() command fails. One fix is to explicitely clear the ifstream, or 00343 // scope the variables. 00344 dataFile.close(); 00345 dataFile.clear(); 00346 00347 00348 //D2: flat array with the values of D2 applied to basis functions. Multi-index is (F,P,D2cardinality) 00349 std::vector<double> basisD2; 00350 fileName = "./testdata/HEX_I2_D2Vals.dat"; 00351 dataFile.open(fileName.c_str()); 00352 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error, 00353 ">>> ERROR (HGRAD_HEX_I2/test01): could not open D2 values data file, test aborted."); 00354 while (!dataFile.eof() ){ 00355 double temp; 00356 string line; // string for one line of input file 00357 std::getline(dataFile, line); // get next line from file 00358 stringstream data_line(line); // convert to stringstream 00359 while(data_line >> temp){ // extract value from line 00360 basisD2.push_back(temp); // push into vector 00361 } 00362 } 00363 dataFile.close(); 00364 dataFile.clear(); 00365 00366 00367 //D3: flat array with the values of D3 applied to basis functions. Multi-index is (F,P,D3cardinality) 00368 std::vector<double> basisD3; 00369 00370 fileName = "./testdata/HEX_I2_D3Vals.dat"; 00371 dataFile.open(fileName.c_str()); 00372 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error, 00373 ">>> ERROR (HGRAD_HEX_I2/test01): could not open D3 values data file, test aborted."); 00374 00375 while (!dataFile.eof() ){ 00376 double temp; 00377 string line; // string for one line of input file 00378 std::getline(dataFile, line); // get next line from file 00379 stringstream data_line(line); // convert to stringstream 00380 while(data_line >> temp){ // extract value from line 00381 basisD3.push_back(temp); // push into vector 00382 } 00383 } 00384 dataFile.close(); 00385 dataFile.clear(); 00386 00387 00388 00389 try{ 00390 00391 // Dimensions for the output arrays: 00392 int numFields = hexBasis.getCardinality(); 00393 int numPoints = hexNodes.dimension(0); 00394 int spaceDim = hexBasis.getBaseCellTopology().getDimension(); 00395 00396 // Generic array for values, grads, curls, etc. that will be properly sized before each call 00397 FieldContainer<double> vals; 00398 00399 // Check VALUE of basis functions: resize vals to rank-2 container: 00400 vals.resize(numFields, numPoints); 00401 hexBasis.getValues(vals, hexNodes, OPERATOR_VALUE); 00402 for (int i = 0; i < numFields; i++) { 00403 for (int j = 0; j < numPoints; j++) { 00404 int l = i + j * numFields; 00405 if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) { 00406 errorFlag++; 00407 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00408 00409 // Output the multi-index of the value where the error is: 00410 *outStream << " At multi-index { "; 00411 *outStream << i << " ";*outStream << j << " "; 00412 *outStream << "} computed value: " << vals(i,j) 00413 << " but reference value: " << basisValues[l] << "\n"; 00414 } 00415 } 00416 } 00417 00418 00419 // Check GRAD of basis function: resize vals to rank-3 container 00420 vals.resize(numFields, numPoints, spaceDim); 00421 hexBasis.getValues(vals, hexNodes, OPERATOR_GRAD); 00422 for (int i = 0; i < numFields; i++) { 00423 for (int j = 0; j < numPoints; j++) { 00424 for (int k = 0; k < spaceDim; k++) { 00425 00426 // basisGrads is (F,P,D), compute offset: 00427 int l = k + j * spaceDim + i * spaceDim * numPoints; 00428 00429 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00430 errorFlag++; 00431 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00432 00433 // Output the multi-index of the value where the error is: 00434 *outStream << " At multi-index { "; 00435 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00436 *outStream << "} computed grad component: " << vals(i,j,k) 00437 << " but reference grad component: " << basisGrads[l] << "\n"; 00438 } 00439 00440 } 00441 } 00442 } 00443 00444 // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD) 00445 hexBasis.getValues(vals, hexNodes, OPERATOR_D1); 00446 for (int i = 0; i < numFields; i++) { 00447 for (int j = 0; j < numPoints; j++) { 00448 for (int k = 0; k < spaceDim; k++) { 00449 00450 // basisGrads is (F,P,D), compute offset: 00451 int l = k + j * spaceDim + i * spaceDim * numPoints; 00452 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00453 errorFlag++; 00454 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00455 00456 // Output the multi-index of the value where the error is: 00457 *outStream << " At multi-index { "; 00458 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00459 *outStream << "} computed D1 component: " << vals(i,j,k) 00460 << " but reference D1 component: " << basisGrads[l] << "\n"; 00461 } 00462 00463 } 00464 } 00465 } 00466 00467 00468 // Check D2 of basis function 00469 int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim); 00470 vals.resize(numFields, numPoints, D2cardinality); 00471 hexBasis.getValues(vals, hexNodes, OPERATOR_D2); 00472 for (int i = 0; i < numFields; i++) { 00473 for (int j = 0; j < numPoints; j++) { 00474 for (int k = 0; k < D2cardinality; k++) { 00475 00476 // basisD2 is (F,P,Dk), compute offset: 00477 int l = k + j * D2cardinality + i * D2cardinality * numPoints; 00478 if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) { 00479 errorFlag++; 00480 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00481 00482 // Output the multi-index of the value where the error is: 00483 *outStream << " At multi-index { "; 00484 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00485 *outStream << "} computed D2 component: " << vals(i,j,k) 00486 << " but reference D2 component: " << basisD2[l] << "\n"; 00487 } 00488 } 00489 } 00490 } 00491 00492 00493 // Check D3 of basis function 00494 int D3cardinality = Intrepid::getDkCardinality(OPERATOR_D3, spaceDim); 00495 vals.resize(numFields, numPoints, D3cardinality); 00496 hexBasis.getValues(vals, hexNodes, OPERATOR_D3); 00497 00498 for (int i = 0; i < numFields; i++) { 00499 for (int j = 0; j < numPoints; j++) { 00500 for (int k = 0; k < D3cardinality; k++) { 00501 00502 // basisD3 is (F,P,Dk), compute offset: 00503 int l = k + j * D3cardinality + i * D3cardinality * numPoints; 00504 if (std::abs(vals(i,j,k) - basisD3[l]) > INTREPID_TOL) { 00505 errorFlag++; 00506 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00507 00508 // Output the multi-index of the value where the error is: 00509 *outStream << " At multi-index { "; 00510 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00511 *outStream << "} computed D3 component: " << vals(i,j,k) 00512 << " but reference D3 component: " << basisD3[l] << "\n"; 00513 } 00514 } 00515 } 00516 } 00517 00518 } 00519 00520 // Catch unexpected errors 00521 catch (std::logic_error err) { 00522 *outStream << err.what() << "\n\n"; 00523 errorFlag = -1000; 00524 }; 00525 00526 *outStream \ 00527 << "\n" 00528 << "===============================================================================\n"\ 00529 << "| TEST 4: correctness of DoF locations |\n"\ 00530 << "===============================================================================\n"; 00531 00532 try{ 00533 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis = 00534 Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM<double, FieldContainer<double> >); 00535 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface = 00536 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis); 00537 00538 FieldContainer<double> cvals; 00539 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality()); 00540 00541 // Check exceptions. 00542 #ifdef HAVE_INTREPID_DEBUG 00543 cvals.resize(1,2,3); 00544 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00545 cvals.resize(3,2); 00546 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00547 cvals.resize(20,2); 00548 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00549 #endif 00550 cvals.resize(20,3); 00551 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--; 00552 // Check if number of thrown exceptions matches the one we expect 00553 if (throwCounter != nException) { 00554 errorFlag++; 00555 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00556 } 00557 00558 // Check mathematical correctness. 00559 basis->getValues(bvals, cvals, OPERATOR_VALUE); 00560 char buffer[120]; 00561 for (int i=0; i<bvals.dimension(0); i++) { 00562 for (int j=0; j<bvals.dimension(1); j++) { 00563 if ((i != j) && (std::abs(bvals(i,j) - 0.0) > INTREPID_TOL)) { 00564 errorFlag++; 00565 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), bvals(i,j), 0.0); 00566 *outStream << buffer; 00567 } 00568 else if ((i == j) && (std::abs(bvals(i,j) - 1.0) > INTREPID_TOL)) { 00569 errorFlag++; 00570 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), bvals(i,j), 1.0); 00571 *outStream << buffer; 00572 } 00573 } 00574 } 00575 00576 } 00577 catch (std::logic_error err){ 00578 *outStream << err.what() << "\n\n"; 00579 errorFlag = -1000; 00580 }; 00581 00582 if (errorFlag != 0) 00583 std::cout << "End Result: TEST FAILED\n"; 00584 else 00585 std::cout << "End Result: TEST PASSED\n"; 00586 00587 // reset format state of std::cout 00588 std::cout.copyfmt(oldFormatState); 00589 00590 return errorFlag; 00591 }
1.7.6.1