|
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_TET_C2_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_TET_C2_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_TET_C2_FEM<double, FieldContainer<double> > tetBasis; 00110 int errorFlag = 0; 00111 00112 // Initialize throw counter for exception testing 00113 int nException = 0; 00114 int throwCounter = 0; 00115 00116 // Define array containing the 10 nodes of Tetrahedron<10>: 4 vertices + 6 edge midpoints 00117 FieldContainer<double> tetNodes(10, 3); 00118 tetNodes(0,0) = 0.0; tetNodes(0,1) = 0.0; tetNodes(0,2) = 0.0; 00119 tetNodes(1,0) = 1.0; tetNodes(1,1) = 0.0; tetNodes(1,2) = 0.0; 00120 tetNodes(2,0) = 0.0; tetNodes(2,1) = 1.0; tetNodes(2,2) = 0.0; 00121 tetNodes(3,0) = 0.0; tetNodes(3,1) = 0.0; tetNodes(3,2) = 1.0; 00122 00123 tetNodes(4,0) = 0.5; tetNodes(4,1) = 0.0; tetNodes(4,2) = 0.0; 00124 tetNodes(5,0) = 0.5; tetNodes(5,1) = 0.5; tetNodes(5,2) = 0.0; 00125 tetNodes(6,0) = 0.0; tetNodes(6,1) = 0.5; tetNodes(6,2) = 0.0; 00126 tetNodes(7,0) = 0.0; tetNodes(7,1) = 0.0; tetNodes(7,2) = 0.5; 00127 tetNodes(8,0) = 0.5; tetNodes(8,1) = 0.0; tetNodes(8,2) = 0.5; 00128 tetNodes(9,0) = 0.0; tetNodes(9,1) = 0.5; tetNodes(9,2) = 0.5; 00129 00130 00131 // Generic array for the output values; needs to be properly resized depending on the operator type 00132 FieldContainer<double> vals; 00133 00134 try{ 00135 // exception #1: CURL cannot be applied to scalar functions 00136 // resize vals to rank-3 container with dimensions (num. points, num. basis functions, arbitrary) 00137 vals.resize(tetBasis.getCardinality(), tetNodes.dimension(0), 4 ); 00138 INTREPID_TEST_COMMAND( tetBasis.getValues(vals, tetNodes, OPERATOR_CURL), throwCounter, nException ); 00139 00140 // exception #2: DIV cannot be applied to scalar functions 00141 // resize vals to rank-2 container with dimensions (num. points, num. basis functions) 00142 vals.resize(tetBasis.getCardinality(), tetNodes.dimension(0) ); 00143 INTREPID_TEST_COMMAND( tetBasis.getValues(vals, tetNodes, OPERATOR_DIV), throwCounter, nException ); 00144 00145 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00146 // getDofTag() to access invalid array elements thereby causing bounds check exception 00147 // exception #3 00148 INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00149 // exception #4 00150 INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00151 // exception #5 00152 INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(0,4,0), throwCounter, nException ); 00153 // exception #6 00154 INTREPID_TEST_COMMAND( tetBasis.getDofTag(10), throwCounter, nException ); 00155 // exception #7 00156 INTREPID_TEST_COMMAND( tetBasis.getDofTag(-1), throwCounter, nException ); 00157 00158 #ifdef HAVE_INTREPID_DEBUG 00159 // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays 00160 // exception #8: input points array must be of rank-2 00161 FieldContainer<double> badPoints1(4, 5, 3); 00162 INTREPID_TEST_COMMAND( tetBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00163 00164 // exception #9 dimension 1 in the input point array must equal space dimension of the cell 00165 FieldContainer<double> badPoints2(4, tetBasis.getBaseCellTopology().getDimension() - 1); 00166 INTREPID_TEST_COMMAND( tetBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00167 00168 // exception #10 output values must be of rank-2 for OPERATOR_VALUE 00169 FieldContainer<double> badVals1(4, 3, 1); 00170 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals1, tetNodes, OPERATOR_VALUE), throwCounter, nException ); 00171 00172 // exception #11 output values must be of rank-3 for OPERATOR_GRAD 00173 FieldContainer<double> badVals2(4, 3); 00174 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals2, tetNodes, OPERATOR_GRAD), throwCounter, nException ); 00175 00176 // exception #12 output values must be of rank-3 for OPERATOR_D1 00177 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals2, tetNodes, OPERATOR_D1), throwCounter, nException ); 00178 00179 // exception #13 output values must be of rank-3 for OPERATOR_D2 00180 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals2, tetNodes, OPERATOR_D2), throwCounter, nException ); 00181 00182 // exception #14 incorrect 0th dimension of output array (must equal number of basis functions) 00183 FieldContainer<double> badVals3(tetBasis.getCardinality() + 1, tetNodes.dimension(0)); 00184 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals3, tetNodes, OPERATOR_VALUE), throwCounter, nException ); 00185 00186 // exception #15 incorrect 1st dimension of output array (must equal number of points) 00187 FieldContainer<double> badVals4(tetBasis.getCardinality(), tetNodes.dimension(0) + 1); 00188 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals4, tetNodes, OPERATOR_VALUE), throwCounter, nException ); 00189 00190 // exception #16: incorrect 2nd dimension of output array (must equal the space dimension) 00191 FieldContainer<double> badVals5(tetBasis.getCardinality(), tetNodes.dimension(0), tetBasis.getBaseCellTopology().getDimension() + 1); 00192 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals5, tetNodes, OPERATOR_GRAD), throwCounter, nException ); 00193 00194 // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 2D) 00195 FieldContainer<double> badVals6(tetBasis.getCardinality(), tetNodes.dimension(0), 40); 00196 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals6, tetNodes, OPERATOR_D1), throwCounter, nException ); 00197 00198 // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 2D) 00199 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals6, tetNodes, OPERATOR_D2), throwCounter, nException ); 00200 #endif 00201 00202 } 00203 catch (std::logic_error err) { 00204 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00205 *outStream << err.what() << '\n'; 00206 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00207 errorFlag = -1000; 00208 }; 00209 00210 // Check if number of thrown exceptions matches the one we expect 00211 if (throwCounter != nException) { 00212 errorFlag++; 00213 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00214 } 00215 00216 *outStream \ 00217 << "\n" 00218 << "===============================================================================\n"\ 00219 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00220 << "===============================================================================\n"; 00221 00222 try{ 00223 std::vector<std::vector<int> > allTags = tetBasis.getAllDofTags(); 00224 00225 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00226 for (unsigned i = 0; i < allTags.size(); i++) { 00227 int bfOrd = tetBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00228 00229 std::vector<int> myTag = tetBasis.getDofTag(bfOrd); 00230 if( !( (myTag[0] == allTags[i][0]) && 00231 (myTag[1] == allTags[i][1]) && 00232 (myTag[2] == allTags[i][2]) && 00233 (myTag[3] == allTags[i][3]) ) ) { 00234 errorFlag++; 00235 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00236 *outStream << " getDofOrdinal( {" 00237 << allTags[i][0] << ", " 00238 << allTags[i][1] << ", " 00239 << allTags[i][2] << ", " 00240 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00241 *outStream << " getDofTag(" << bfOrd << ") = { " 00242 << myTag[0] << ", " 00243 << myTag[1] << ", " 00244 << myTag[2] << ", " 00245 << myTag[3] << "}\n"; 00246 } 00247 } 00248 00249 // Now do the same but loop over basis functions 00250 for( int bfOrd = 0; bfOrd < tetBasis.getCardinality(); bfOrd++) { 00251 std::vector<int> myTag = tetBasis.getDofTag(bfOrd); 00252 int myBfOrd = tetBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00253 if( bfOrd != myBfOrd) { 00254 errorFlag++; 00255 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00256 *outStream << " getDofTag(" << bfOrd << ") = { " 00257 << myTag[0] << ", " 00258 << myTag[1] << ", " 00259 << myTag[2] << ", " 00260 << myTag[3] << "} but getDofOrdinal({" 00261 << myTag[0] << ", " 00262 << myTag[1] << ", " 00263 << myTag[2] << ", " 00264 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00265 } 00266 } 00267 } 00268 catch (std::logic_error err){ 00269 *outStream << err.what() << "\n\n"; 00270 errorFlag = -1000; 00271 }; 00272 00273 *outStream \ 00274 << "\n" 00275 << "===============================================================================\n"\ 00276 << "| TEST 3: correctness of basis function values |\n"\ 00277 << "===============================================================================\n"; 00278 00279 outStream -> precision(20); 00280 00281 // VALUE: in (F,P) format 00282 double basisValues[] = { 00283 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, \ 00284 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, \ 00285 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, \ 00286 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00287 1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 0, 0, 0, 0, 0, 0, 0, \ 00288 0, 0, 0, 1.00000 }; 00289 00290 // GRAD and D1: in (F,P,D) format 00291 double basisGrads[] = { 00292 -3.00000, -3.00000, -3.00000, 1.00000, 1.00000, 1.00000, 1.00000, \ 00293 1.00000, 1.00000, 1.00000, 1.00000, 1.00000, -1.00000, -1.00000, \ 00294 -1.00000, 1.00000, 1.00000, 1.00000, -1.00000, -1.00000, -1.00000, \ 00295 -1.00000, -1.00000, -1.00000, 1.00000, 1.00000, 1.00000, 1.00000, \ 00296 1.00000, 1.00000, -1.00000, 0, 0, 3.00000, 0, 0, -1.00000, 0, 0, \ 00297 -1.00000, 0, 0, 1.00000, 0, 0, 1.00000, 0, 0, -1.00000, 0, 0, \ 00298 -1.00000, 0, 0, 1.00000, 0, 0, -1.00000, 0, 0, 0, -1.00000, 0, 0, \ 00299 -1.00000, 0, 0, 3.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \ 00300 1.00000, 0, 0, 1.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \ 00301 1.00000, 0, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \ 00302 3.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, -1.00000, 0, 0, \ 00303 1.00000, 0, 0, 1.00000, 0, 0, 1.00000, 4.00000, 0, 0, -4.00000, \ 00304 -4.00000, -4.00000, 0, 0, 0, 0, 0, 0, 0, -2.00000, -2.00000, \ 00305 -2.00000, -2.00000, -2.00000, 2.00000, 0, 0, 2.00000, 0, 0, -2.00000, \ 00306 -2.00000, -2.00000, 0, 0, 0, 0, 0, 0, 0, 4.00000, 0, 4.00000, 0, 0, \ 00307 0, 0, 0, 0, 2.00000, 0, 2.00000, 2.00000, 0, 2.00000, 0, 0, 0, 0, 0, \ 00308 0, 2.00000, 0, 2.00000, 0, 0, 0, 4.00000, 0, 0, 0, 0, -4.00000, \ 00309 -4.00000, -4.00000, 0, 0, 0, 0, 2.00000, 0, -2.00000, -2.00000, \ 00310 -2.00000, -2.00000, 0, -2.00000, 0, 2.00000, 0, 0, 0, 0, -2.00000, \ 00311 -2.00000, -2.00000, 0, 0, 4.00000, 0, 0, 0, 0, 0, 0, -4.00000, \ 00312 -4.00000, -4.00000, 0, 0, 2.00000, 0, 0, 0, 0, 0, 2.00000, -2.00000, \ 00313 -2.00000, 0, -2.00000, -2.00000, -2.00000, -2.00000, -2.00000, \ 00314 -2.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 4.00000, 0, 0, 0, 0, \ 00315 2.00000, 0, 0, 2.00000, 0, 0, 0, 2.00000, 0, 0, 2.00000, 0, 2.00000, \ 00316 2.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 4.00000, 0, 4.00000, 0, 0, 0, \ 00317 0, 0, 0, 2.00000, 0, 0, 2.00000, 0, 2.00000, 0, 0, 2.00000, 0, 0, \ 00318 2.00000, 2.00000}; 00319 00320 // D2 values in (F,P, Dk) format 00321 double basisD2[]={ 00322 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00323 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00324 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00325 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00326 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00327 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00328 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00329 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, \ 00330 4.00000, 4.00000, 4.00000, 4.00000, 4.00000, 0, 0, 0, 0, 0, 4.00000, \ 00331 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, \ 00332 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \ 00333 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \ 00334 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, \ 00335 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \ 00336 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \ 00337 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 0, 0, 4.00000, \ 00338 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, \ 00339 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \ 00340 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \ 00341 0, 0, 4.00000, -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, \ 00342 -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, \ 00343 -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, \ 00344 0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, \ 00345 -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, -8.00000, \ 00346 -4.00000, -4.00000, 0, 0, 0, -8.00000, -4.00000, -4.00000, 0, 0, 0, \ 00347 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \ 00348 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, \ 00349 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, \ 00350 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, -4.00000, 0, -8.00000, -4.00000, \ 00351 0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, -8.00000, \ 00352 -4.00000, 0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, \ 00353 -8.00000, -4.00000, 0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, \ 00354 -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, -8.00000, \ 00355 -4.00000, 0, 0, -4.00000, 0, -8.00000, -4.00000, 0, 0, -4.00000, 0, \ 00356 -8.00000, -4.00000, 0, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, \ 00357 -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, -4.00000, \ 00358 -8.00000, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, \ 00359 -4.00000, -8.00000, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, \ 00360 -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, -4.00000, \ 00361 -8.00000, 0, 0, -4.00000, 0, -4.00000, -8.00000, 0, 0, -4.00000, 0, \ 00362 -4.00000, -8.00000, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \ 00363 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, \ 00364 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, \ 00365 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 0, 0, \ 00366 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, \ 00367 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, \ 00368 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, 0, 0, 0, 4.00000, 0, 0, \ 00369 0, 0, 0, 4.00000, 0 00370 }; 00371 00372 00373 try{ 00374 00375 // Dimensions for the output arrays: 00376 int numFields = tetBasis.getCardinality(); 00377 int numPoints = tetNodes.dimension(0); 00378 int spaceDim = tetBasis.getBaseCellTopology().getDimension(); 00379 00380 // Generic array for values, grads, curls, etc. that will be properly sized before each call 00381 FieldContainer<double> vals; 00382 00383 // Check VALUE of basis functions: resize vals to rank-2 container: 00384 vals.resize(numFields, numPoints); 00385 tetBasis.getValues(vals, tetNodes, OPERATOR_VALUE); 00386 for (int i = 0; i < numFields; i++) { 00387 for (int j = 0; j < numPoints; j++) { 00388 int l = i + j * numFields; 00389 if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) { 00390 errorFlag++; 00391 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00392 00393 // Output the multi-index of the value where the error is: 00394 *outStream << " At multi-index { "; 00395 *outStream << i << " ";*outStream << j << " "; 00396 *outStream << "} computed value: " << vals(i,j) 00397 << " but reference value: " << basisValues[l] << "\n"; 00398 } 00399 } 00400 } 00401 00402 // Check GRAD of basis function: resize vals to rank-3 container 00403 vals.resize(numFields, numPoints, spaceDim); 00404 tetBasis.getValues(vals, tetNodes, OPERATOR_GRAD); 00405 for (int i = 0; i < numFields; i++) { 00406 for (int j = 0; j < numPoints; j++) { 00407 for (int k = 0; k < spaceDim; k++) { 00408 00409 // basisGrads is (F,P,D), compute offset: 00410 int l = k + j * spaceDim + i * spaceDim * numPoints; 00411 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00412 errorFlag++; 00413 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00414 00415 // Output the multi-index of the value where the error is: 00416 *outStream << " At multi-index { "; 00417 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00418 *outStream << "} computed grad component: " << vals(i,j,k) 00419 << " but reference grad component: " << basisGrads[l] << "\n"; 00420 } 00421 } 00422 } 00423 } 00424 00425 // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD) 00426 tetBasis.getValues(vals, tetNodes, OPERATOR_D1); 00427 for (int i = 0; i < numFields; i++) { 00428 for (int j = 0; j < numPoints; j++) { 00429 for (int k = 0; k < spaceDim; k++) { 00430 00431 // basisGrads is (F,P,D), compute offset: 00432 int l = k + j * spaceDim + i * spaceDim * numPoints; 00433 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00434 errorFlag++; 00435 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00436 00437 // Output the multi-index of the value where the error is: 00438 *outStream << " At multi-index { "; 00439 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00440 *outStream << "} computed D1 component: " << vals(i,j,k) 00441 << " but reference D1 component: " << basisGrads[l] << "\n"; 00442 } 00443 } 00444 } 00445 } 00446 00447 // Check D2 of basis function 00448 int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim); 00449 vals.resize(numFields, numPoints, D2cardinality); 00450 tetBasis.getValues(vals, tetNodes, OPERATOR_D2); 00451 for (int i = 0; i < numFields; i++) { 00452 for (int j = 0; j < numPoints; j++) { 00453 for (int k = 0; k < D2cardinality; k++) { 00454 00455 // basisD2 is (F,P,Dk), compute offset: 00456 int l = k + j * D2cardinality + i * D2cardinality * numPoints; 00457 if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) { 00458 errorFlag++; 00459 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00460 00461 // Output the multi-index of the value where the error is: 00462 *outStream << " At multi-index { "; 00463 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00464 *outStream << "} computed D2 component: " << vals(i,j,k) 00465 << " but reference D2 component: " << basisD2[l] << "\n"; 00466 } 00467 } 00468 } 00469 } 00470 00471 00472 // Check all higher derivatives - must be zero. 00473 for(EOperator op = OPERATOR_D3; op < OPERATOR_MAX; op++) { 00474 00475 // The last dimension is the number of kth derivatives and needs to be resized for every Dk 00476 int DkCardin = Intrepid::getDkCardinality(op, spaceDim); 00477 vals.resize(numFields, numPoints, DkCardin); 00478 00479 tetBasis.getValues(vals, tetNodes, op); 00480 for (int i = 0; i < vals.size(); i++) { 00481 if (std::abs(vals[i]) > INTREPID_TOL) { 00482 errorFlag++; 00483 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00484 00485 // Get the multi-index of the value where the error is and the operator order 00486 std::vector<int> myIndex; 00487 vals.getMultiIndex(myIndex,i); 00488 int ord = Intrepid::getOperatorOrder(op); 00489 *outStream << " At multi-index { "; 00490 for(int j = 0; j < vals.rank(); j++) { 00491 *outStream << myIndex[j] << " "; 00492 } 00493 *outStream << "} computed D"<< ord <<" component: " << vals[i] 00494 << " but reference D" << ord << " component: 0 \n"; 00495 } 00496 } 00497 } 00498 } 00499 00500 // Catch unexpected errors 00501 catch (std::logic_error err) { 00502 *outStream << err.what() << "\n\n"; 00503 errorFlag = -1000; 00504 }; 00505 00506 if (errorFlag != 0) 00507 std::cout << "End Result: TEST FAILED\n"; 00508 else 00509 std::cout << "End Result: TEST PASSED\n"; 00510 00511 // reset format state of std::cout 00512 std::cout.copyfmt(oldFormatState); 00513 00514 return errorFlag; 00515 }
1.7.6.1