|
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_HCURL_WEDGE_I1_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_HCURL_WEDGE_I1_FEM) |\n" \ 00093 << "| |\n" \ 00094 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \ 00095 << "| 2) Basis values for VALUE and CURL 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_HCURL_WEDGE_I1_FEM<double, FieldContainer<double> > wedgeBasis; 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 6 vertices of the reference WEDGE and 6 other points. 00117 FieldContainer<double> wedgeNodes(12, 3); 00118 wedgeNodes(0,0) = 0.0; wedgeNodes(0,1) = 0.0; wedgeNodes(0,2) = -1.0; 00119 wedgeNodes(1,0) = 1.0; wedgeNodes(1,1) = 0.0; wedgeNodes(1,2) = -1.0; 00120 wedgeNodes(2,0) = 0.0; wedgeNodes(2,1) = 1.0; wedgeNodes(2,2) = -1.0; 00121 wedgeNodes(3,0) = 0.0; wedgeNodes(3,1) = 0.0; wedgeNodes(3,2) = 1.0; 00122 wedgeNodes(4,0) = 1.0; wedgeNodes(4,1) = 0.0; wedgeNodes(4,2) = 1.0; 00123 wedgeNodes(5,0) = 0.0; wedgeNodes(5,1) = 1.0; wedgeNodes(5,2) = 1.0; 00124 00125 wedgeNodes(6,0) = 0.25; wedgeNodes(6,1) = 0.5; wedgeNodes(6,2) = -1.0; 00126 wedgeNodes(7,0) = 0.5; wedgeNodes(7,1) = 0.25; wedgeNodes(7,2) = 0.0; 00127 wedgeNodes(8,0) = 0.25; wedgeNodes(8,1) = 0.25; wedgeNodes(8,2) = 1.0; 00128 wedgeNodes(9,0) = 0.25; wedgeNodes(9,1) = 0.0; wedgeNodes(9,2) = 0.75; 00129 wedgeNodes(10,0)= 0.0; wedgeNodes(10,1)= 0.5; wedgeNodes(10,2)= -0.25; 00130 wedgeNodes(11,0)= 0.5; wedgeNodes(11,1)= 0.5; wedgeNodes(11,2)= 0.0; 00131 00132 00133 00134 // Generic array for the output values; needs to be properly resized depending on the operator type 00135 FieldContainer<double> vals; 00136 00137 try{ 00138 // exception #1: GRAD cannot be applied to HCURL functions 00139 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00140 vals.resize(wedgeBasis.getCardinality(), wedgeNodes.dimension(0), 3 ); 00141 INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_GRAD), throwCounter, nException ); 00142 00143 // exception #2: DIV cannot be applied to HCURL functions 00144 // resize vals to rank-2 container with dimensions (num. basis functions, num. points) 00145 vals.resize(wedgeBasis.getCardinality(), wedgeNodes.dimension(0)); 00146 INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_DIV), throwCounter, nException ); 00147 00148 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00149 // getDofTag() to access invalid array elements thereby causing bounds check exception 00150 // exception #3 00151 INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00152 // exception #4 00153 INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00154 // exception #5 00155 INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(0,4,1), throwCounter, nException ); 00156 // exception #6 00157 INTREPID_TEST_COMMAND( wedgeBasis.getDofTag(10), throwCounter, nException ); 00158 // exception #7 00159 INTREPID_TEST_COMMAND( wedgeBasis.getDofTag(-1), throwCounter, nException ); 00160 00161 #ifdef HAVE_INTREPID_DEBUG 00162 // Exceptions 8- test exception handling with incorrectly dimensioned input/output arrays 00163 // exception #8: input points array must be of rank-2 00164 FieldContainer<double> badPoints1(4, 5, 3); 00165 INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00166 00167 // exception #9 dimension 1 in the input point array must equal space dimension of the cell 00168 FieldContainer<double> badPoints2(4, 2); 00169 INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00170 00171 // exception #10 output values must be of rank-3 for OPERATOR_VALUE 00172 FieldContainer<double> badVals1(4, 3); 00173 INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals1, wedgeNodes, OPERATOR_VALUE), throwCounter, nException ); 00174 00175 // exception #11 output values must be of rank-3 for OPERATOR_CURL 00176 INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals1, wedgeNodes, OPERATOR_CURL), throwCounter, nException ); 00177 00178 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions) 00179 FieldContainer<double> badVals2(wedgeBasis.getCardinality() + 1, wedgeNodes.dimension(0), 3); 00180 INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals2, wedgeNodes, OPERATOR_VALUE), throwCounter, nException ); 00181 00182 // exception #13 incorrect 1st dimension of output array (must equal number of points) 00183 FieldContainer<double> badVals3(wedgeBasis.getCardinality(), wedgeNodes.dimension(0) + 1, 3); 00184 INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals3, wedgeNodes, OPERATOR_VALUE), throwCounter, nException ); 00185 00186 // exception #14: incorrect 2nd dimension of output array (must equal the space dimension) 00187 FieldContainer<double> badVals4(wedgeBasis.getCardinality(), wedgeNodes.dimension(0), 4); 00188 INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals4, wedgeNodes, OPERATOR_VALUE), throwCounter, nException ); 00189 00190 // exception #15: incorrect 2nd dimension of output array (must equal the space dimension) 00191 INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals4, wedgeNodes, OPERATOR_CURL), throwCounter, nException ); 00192 #endif 00193 00194 } 00195 catch (std::logic_error err) { 00196 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00197 *outStream << err.what() << '\n'; 00198 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00199 errorFlag = -1000; 00200 }; 00201 00202 // Check if number of thrown exceptions matches the one we expect 00203 // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match. 00204 if (throwCounter != nException) { 00205 errorFlag++; 00206 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00207 } 00208 00209 *outStream \ 00210 << "\n" 00211 << "===============================================================================\n"\ 00212 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00213 << "===============================================================================\n"; 00214 00215 try{ 00216 std::vector<std::vector<int> > allTags = wedgeBasis.getAllDofTags(); 00217 00218 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00219 for (unsigned i = 0; i < allTags.size(); i++) { 00220 int bfOrd = wedgeBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00221 00222 std::vector<int> myTag = wedgeBasis.getDofTag(bfOrd); 00223 if( !( (myTag[0] == allTags[i][0]) && 00224 (myTag[1] == allTags[i][1]) && 00225 (myTag[2] == allTags[i][2]) && 00226 (myTag[3] == allTags[i][3]) ) ) { 00227 errorFlag++; 00228 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00229 *outStream << " getDofOrdinal( {" 00230 << allTags[i][0] << ", " 00231 << allTags[i][1] << ", " 00232 << allTags[i][2] << ", " 00233 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00234 *outStream << " getDofTag(" << bfOrd << ") = { " 00235 << myTag[0] << ", " 00236 << myTag[1] << ", " 00237 << myTag[2] << ", " 00238 << myTag[3] << "}\n"; 00239 } 00240 } 00241 00242 // Now do the same but loop over basis functions 00243 for( int bfOrd = 0; bfOrd < wedgeBasis.getCardinality(); bfOrd++) { 00244 std::vector<int> myTag = wedgeBasis.getDofTag(bfOrd); 00245 int myBfOrd = wedgeBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00246 if( bfOrd != myBfOrd) { 00247 errorFlag++; 00248 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00249 *outStream << " getDofTag(" << bfOrd << ") = { " 00250 << myTag[0] << ", " 00251 << myTag[1] << ", " 00252 << myTag[2] << ", " 00253 << myTag[3] << "} but getDofOrdinal({" 00254 << myTag[0] << ", " 00255 << myTag[1] << ", " 00256 << myTag[2] << ", " 00257 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00258 } 00259 } 00260 } 00261 catch (std::logic_error err){ 00262 *outStream << err.what() << "\n\n"; 00263 errorFlag = -1000; 00264 }; 00265 00266 *outStream \ 00267 << "\n" 00268 << "===============================================================================\n"\ 00269 << "| TEST 3: correctness of basis function values |\n"\ 00270 << "===============================================================================\n"; 00271 00272 outStream -> precision(20); 00273 00274 // VALUE: Each row pair gives the 9x3 correct basis set values at an evaluation point: (P,F,D) layout 00275 double basisValues[] = { 00276 1.00000, 0, 0, 0, 0, 0, 0, -1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00277 0, 0.500000, 0, 0, 0, 0, 0, 0, 1.00000, 1.00000, 0, 0, 1.00000, 0, 0, \ 00278 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.500000, 0, 0, 0, 0, \ 00279 0, 0, -1.00000, 0, 0, -1.00000, -1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00280 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.500000, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00281 1.00000, 0, 0, 0, 0, 0, 0, -1.00000, 0, 0, 0, 0.500000, 0, 0, 0, 0, \ 00282 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 1.00000, 0, 0, 1.00000, 0, \ 00283 0, 0, 0, 0, 0, 0, 0, 0, 0.500000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \ 00284 0, 0, 0, -1.00000, 0, 0, -1.00000, -1.00000, 0, 0, 0, 0, 0, 0, 0, 0, \ 00285 0, 0.500000, 0.500000, 0.250000, 0, -0.500000, 0.250000, 0, \ 00286 -0.500000, -0.750000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.125000, \ 00287 0, 0, 0.125000, 0, 0, 0.250000, 0.375000, 0.250000, 0, -0.125000, \ 00288 0.250000, 0, -0.125000, -0.250000, 0, 0.375000, 0.250000, 0, \ 00289 -0.125000, 0.250000, 0, -0.125000, -0.250000, 0, 0, 0, 0.125000, 0, \ 00290 0, 0.250000, 0, 0, 0.125000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.750000, \ 00291 0.250000, 0, -0.250000, 0.250000, 0, -0.250000, -0.750000, 0, 0, 0, \ 00292 0.250000, 0, 0, 0.125000, 0, 0, 0.125000, 0.125000, 0.0312500, 0, 0, \ 00293 0.0312500, 0, 0, -0.0937500, 0, 0.875000, 0.218750, 0, 0, 0.218750, \ 00294 0, 0, -0.656250, 0, 0, 0, 0.375000, 0, 0, 0.125000, 0, 0, 0, \ 00295 0.312500, 0, 0, -0.312500, 0, 0, -0.312500, -0.625000, 0, 0.187500, \ 00296 0, 0, -0.187500, 0, 0, -0.187500, -0.375000, 0, 0, 0, 0.250000, 0, 0, \ 00297 0, 0, 0, 0.250000, 0.250000, 0.250000, 0, -0.250000, 0.250000, 0, \ 00298 -0.250000, -0.250000, 0, 0.250000, 0.250000, 0, -0.250000, 0.250000, \ 00299 0, -0.250000, -0.250000, 0, 0, 0, 0, 0, 0, 0.250000, 0, 0, 0.250000 00300 }; 00301 00302 // CURL: each row pair gives the 9x3 correct values of the curls of the 9 basis functions: (P,F,D) layout 00303 double basisCurls[] = { 00304 0, -0.500000, 2.00000, 0, 0, 2.00000, -0.500000, 0, 2.00000, 0, \ 00305 0.500000, 0, 0, 0, 0, 0.500000, 0, 0, -0.500000, 0.500000, 0, 0, \ 00306 -0.500000, 0, 0.500000, 0, 0, 0.500000, -0.500000, 2.00000, 0.500000, \ 00307 0, 2.00000, 0, 0, 2.00000, -0.500000, 0.500000, 0, -0.500000, 0, 0, \ 00308 0, 0, 0, -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0, \ 00309 0, 2.00000, 0, 0.500000, 2.00000, -0.500000, 0.500000, 2.00000, 0, 0, \ 00310 0, 0, -0.500000, 0, 0.500000, -0.500000, 0, -0.500000, 0.500000, 0, \ 00311 0, -0.500000, 0, 0.500000, 0, 0, 0, -0.500000, 0, 0, 0, 0, -0.500000, \ 00312 0, 0, 0, 0.500000, 2.00000, 0, 0, 2.00000, 0.500000, 0, 2.00000, \ 00313 -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0.500000, \ 00314 -0.500000, 0, 0.500000, 0, 0, 0, 0, 0, -0.500000, 0.500000, 2.00000, \ 00315 -0.500000, 0, 2.00000, 0, 0, 2.00000, -0.500000, 0.500000, 0, 0, \ 00316 -0.500000, 0, 0.500000, 0, 0, 0, 0, 0, 0, 0.500000, 0, -0.500000, \ 00317 0.500000, 0, 0, 0, 2.00000, 0, -0.500000, 2.00000, 0.500000, \ 00318 -0.500000, 2.00000, -0.500000, 0.500000, 0, 0, -0.500000, 0, \ 00319 0.500000, 0, 0, 0.125000, -0.250000, 2.00000, 0.125000, 0.250000, \ 00320 2.00000, -0.375000, 0.250000, 2.00000, -0.125000, 0.250000, 0, \ 00321 -0.125000, -0.250000, 0, 0.375000, -0.250000, 0, -0.500000, 0.500000, \ 00322 0, 0, -0.500000, 0, 0.500000, 0, 0, 0.250000, -0.375000, 1.00000, \ 00323 0.250000, 0.125000, 1.00000, -0.250000, 0.125000, 1.00000, -0.250000, \ 00324 0.375000, 1.00000, -0.250000, -0.125000, 1.00000, 0.250000, \ 00325 -0.125000, 1.00000, -0.500000, 0.500000, 0, 0, -0.500000, 0, \ 00326 0.500000, 0, 0, 0.125000, -0.375000, 0, 0.125000, 0.125000, 0, \ 00327 -0.375000, 0.125000, 0, -0.125000, 0.375000, 2.00000, -0.125000, \ 00328 -0.125000, 2.00000, 0.375000, -0.125000, 2.00000, -0.500000, \ 00329 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0.125000, -0.500000, \ 00330 0.250000, 0.125000, 0, 0.250000, -0.375000, 0, 0.250000, -0.125000, \ 00331 0.500000, 1.75000, -0.125000, 0, 1.75000, 0.375000, 0, 1.75000, \ 00332 -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0, \ 00333 -0.250000, 1.25000, 0, 0.250000, 1.25000, -0.500000, 0.250000, \ 00334 1.25000, 0, 0.250000, 0.750000, 0, -0.250000, 0.750000, 0.500000, \ 00335 -0.250000, 0.750000, -0.500000, 0.500000, 0, 0, -0.500000, 0, \ 00336 0.500000, 0, 0, 0.250000, -0.250000, 1.00000, 0.250000, 0.250000, \ 00337 1.00000, -0.250000, 0.250000, 1.00000, -0.250000, 0.250000, 1.00000, \ 00338 -0.250000, -0.250000, 1.00000, 0.250000, -0.250000, 1.00000, \ 00339 -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0 00340 }; 00341 00342 try{ 00343 00344 // Dimensions for the output arrays: 00345 int numFields = wedgeBasis.getCardinality(); 00346 int numPoints = wedgeNodes.dimension(0); 00347 int spaceDim = wedgeBasis.getBaseCellTopology().getDimension(); 00348 00349 // Generic array for values and curls that will be properly sized before each call 00350 FieldContainer<double> vals; 00351 00352 // Check VALUE of basis functions: resize vals to rank-3 container: 00353 vals.resize(numFields, numPoints, spaceDim); 00354 wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_VALUE); 00355 for (int i = 0; i < numFields; i++) { 00356 for (int j = 0; j < numPoints; j++) { 00357 for (int k = 0; k < spaceDim; k++) { 00358 00359 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k 00360 int l = k + i * spaceDim + j * spaceDim * numFields; 00361 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) { 00362 errorFlag++; 00363 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00364 00365 // Output the multi-index of the value where the error is: 00366 *outStream << " At multi-index { "; 00367 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00368 *outStream << "} computed value: " << vals(i,j,k) 00369 << " but reference value: " << basisValues[l] << "\n"; 00370 } 00371 } 00372 } 00373 } 00374 00375 // Check CURL of basis function: resize vals to rank-3 container 00376 vals.resize(numFields, numPoints, spaceDim); 00377 wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_CURL); 00378 for (int i = 0; i < numFields; i++) { 00379 for (int j = 0; j < numPoints; j++) { 00380 for (int k = 0; k < spaceDim; k++) { 00381 00382 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k 00383 int l = k + i * spaceDim + j * spaceDim * numFields; 00384 if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) { 00385 errorFlag++; 00386 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00387 00388 // Output the multi-index of the value where the error is: 00389 *outStream << " At multi-index { "; 00390 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00391 *outStream << "} computed curl component: " << vals(i,j,k) 00392 << " but reference curl component: " << basisCurls[l] << "\n"; 00393 } 00394 } 00395 } 00396 } 00397 00398 } 00399 00400 // Catch unexpected errors 00401 catch (std::logic_error err) { 00402 *outStream << err.what() << "\n\n"; 00403 errorFlag = -1000; 00404 }; 00405 00406 if (errorFlag != 0) 00407 std::cout << "End Result: TEST FAILED\n"; 00408 else 00409 std::cout << "End Result: TEST PASSED\n"; 00410 00411 // reset format state of std::cout 00412 std::cout.copyfmt(oldFormatState); 00413 00414 return errorFlag; 00415 }
1.7.6.1