|
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_TET_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_TET_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) or |\n" \ 00098 << "| Denis Ridzal (dridzal@sandia.gov). |\n" \ 00099 << "| |\n" \ 00100 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00101 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00102 << "| |\n" \ 00103 << "===============================================================================\n"\ 00104 << "| TEST 1: Basis creation, exception testing |\n"\ 00105 << "===============================================================================\n"; 00106 00107 // Define basis and error flag 00108 Basis_HCURL_TET_I1_FEM<double, FieldContainer<double> > tetBasis; 00109 int errorFlag = 0; 00110 00111 // Define throw number for exception testing 00112 int nException = 0; 00113 int throwCounter = 0; 00114 00115 // Define array containing the 4 vertices of the reference TET and its 6 edge centers. 00116 FieldContainer<double> tetNodes(10, 3); 00117 tetNodes(0,0) = 0.0; tetNodes(0,1) = 0.0; tetNodes(0,2) = 0.0; 00118 tetNodes(1,0) = 1.0; tetNodes(1,1) = 0.0; tetNodes(1,2) = 0.0; 00119 tetNodes(2,0) = 0.0; tetNodes(2,1) = 1.0; tetNodes(2,2) = 0.0; 00120 tetNodes(3,0) = 0.0; tetNodes(3,1) = 0.0; tetNodes(3,2) = 1.0; 00121 tetNodes(4,0) = 0.5; tetNodes(4,1) = 0.0; tetNodes(4,2) = 0.0; 00122 tetNodes(5,0) = 0.5; tetNodes(5,1) = 0.5; tetNodes(5,2) = 0.0; 00123 tetNodes(6,0) = 0.0; tetNodes(6,1) = 0.5; tetNodes(6,2) = 0.0; 00124 tetNodes(7,0) = 0.0; tetNodes(7,1) = 0.0; tetNodes(7,2) = 0.5; 00125 tetNodes(8,0) = 0.5; tetNodes(8,1) = 0.0; tetNodes(8,2) = 0.5; 00126 tetNodes(9,0) = 0.0; tetNodes(9,1) = 0.5; tetNodes(9,2) = 0.5; 00127 00128 00129 // Generic array for the output values; needs to be properly resized depending on the operator type 00130 FieldContainer<double> vals; 00131 00132 00133 try{ 00134 // exception #1: GRAD cannot be applied to HCURL functions 00135 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00136 vals.resize(tetBasis.getCardinality(), tetNodes.dimension(0), 3 ); 00137 INTREPID_TEST_COMMAND( tetBasis.getValues(vals, tetNodes, OPERATOR_GRAD), throwCounter, nException ); 00138 00139 // exception #2: DIV cannot be applied to HCURL functions 00140 // resize vals to rank-2 container with dimensions (num. basis functions, num. points) 00141 vals.resize(tetBasis.getCardinality(), tetNodes.dimension(0) ); 00142 INTREPID_TEST_COMMAND( tetBasis.getValues(vals, tetNodes, OPERATOR_DIV), throwCounter, nException ); 00143 00144 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00145 // getDofTag() to access invalid array elements thereby causing bounds check exception 00146 // exception #3 00147 INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00148 // exception #4 00149 INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00150 // exception #5 00151 INTREPID_TEST_COMMAND( tetBasis.getDofOrdinal(0,4,1), throwCounter, nException ); 00152 // exception #6 00153 INTREPID_TEST_COMMAND( tetBasis.getDofTag(7), throwCounter, nException ); 00154 // exception #7 00155 INTREPID_TEST_COMMAND( tetBasis.getDofTag(-1), throwCounter, nException ); 00156 00157 #ifdef HAVE_INTREPID_DEBUG 00158 // Exceptions 8-15 test exception handling with incorrectly dimensioned input/output arrays 00159 // exception #8: input points array must be of rank-2 00160 FieldContainer<double> badPoints1(4, 5, 3); 00161 INTREPID_TEST_COMMAND( tetBasis.getValues(vals,badPoints1,OPERATOR_VALUE), throwCounter, nException ); 00162 00163 // exception #9 dimension 1 in the input point array must equal space dimension of the cell 00164 FieldContainer<double> badPoints2(4, 2); 00165 INTREPID_TEST_COMMAND( tetBasis.getValues(vals,badPoints2,OPERATOR_VALUE), throwCounter, nException ); 00166 00167 // exception #10 output values must be of rank-3 for OPERATOR_VALUE 00168 FieldContainer<double> badVals1(4, 3); 00169 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals1,tetNodes,OPERATOR_VALUE), throwCounter, nException ); 00170 00171 // exception #11 output values must be of rank-3 for OPERATOR_CURL 00172 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals1,tetNodes,OPERATOR_CURL), throwCounter, nException ); 00173 00174 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions) 00175 FieldContainer<double> badVals2(tetBasis.getCardinality() + 1, tetNodes.dimension(0), 3); 00176 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals2,tetNodes,OPERATOR_VALUE), throwCounter, nException ); 00177 00178 // exception #13 incorrect 1st dimension of output array (must equal number of points) 00179 FieldContainer<double> badVals3(tetBasis.getCardinality(), tetNodes.dimension(0) + 1, 3); 00180 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals3,tetNodes,OPERATOR_VALUE), throwCounter, nException ); 00181 00182 // exception #14: incorrect 2nd dimension of output array (must equal the space dimension) 00183 FieldContainer<double> badVals4(tetBasis.getCardinality(), tetNodes.dimension(0), 4); 00184 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals4,tetNodes,OPERATOR_VALUE), throwCounter, nException ); 00185 00186 // exception #15: incorrect 2nd dimension of output array (must equal the space dimension) 00187 INTREPID_TEST_COMMAND( tetBasis.getValues(badVals4,tetNodes,OPERATOR_CURL), throwCounter, nException ); 00188 #endif 00189 00190 } 00191 catch (std::logic_error err) { 00192 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00193 *outStream << err.what() << '\n'; 00194 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00195 errorFlag = -1000; 00196 }; 00197 00198 // Check if number of thrown exceptions matches the one we expect 00199 // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match. 00200 if (throwCounter != nException) { 00201 errorFlag++; 00202 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00203 } 00204 00205 *outStream \ 00206 << "\n" 00207 << "===============================================================================\n"\ 00208 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00209 << "===============================================================================\n"; 00210 00211 try{ 00212 std::vector<std::vector<int> > allTags = tetBasis.getAllDofTags(); 00213 00214 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00215 for (unsigned i = 0; i < allTags.size(); i++) { 00216 int bfOrd = tetBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00217 00218 std::vector<int> myTag = tetBasis.getDofTag(bfOrd); 00219 if( !( (myTag[0] == allTags[i][0]) && 00220 (myTag[1] == allTags[i][1]) && 00221 (myTag[2] == allTags[i][2]) && 00222 (myTag[3] == allTags[i][3]) ) ) { 00223 errorFlag++; 00224 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00225 *outStream << " getDofOrdinal( {" 00226 << allTags[i][0] << ", " 00227 << allTags[i][1] << ", " 00228 << allTags[i][2] << ", " 00229 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00230 *outStream << " getDofTag(" << bfOrd << ") = { " 00231 << myTag[0] << ", " 00232 << myTag[1] << ", " 00233 << myTag[2] << ", " 00234 << myTag[3] << "}\n"; 00235 } 00236 } 00237 00238 // Now do the same but loop over basis functions 00239 for( int bfOrd = 0; bfOrd < tetBasis.getCardinality(); bfOrd++) { 00240 std::vector<int> myTag = tetBasis.getDofTag(bfOrd); 00241 int myBfOrd = tetBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00242 if( bfOrd != myBfOrd) { 00243 errorFlag++; 00244 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00245 *outStream << " getDofTag(" << bfOrd << ") = { " 00246 << myTag[0] << ", " 00247 << myTag[1] << ", " 00248 << myTag[2] << ", " 00249 << myTag[3] << "} but getDofOrdinal({" 00250 << myTag[0] << ", " 00251 << myTag[1] << ", " 00252 << myTag[2] << ", " 00253 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00254 } 00255 } 00256 } 00257 catch (std::logic_error err){ 00258 *outStream << err.what() << "\n\n"; 00259 errorFlag = -1000; 00260 }; 00261 00262 *outStream \ 00263 << "\n" 00264 << "===============================================================================\n"\ 00265 << "| TEST 3: correctness of basis function values |\n"\ 00266 << "===============================================================================\n"; 00267 00268 outStream -> precision(20); 00269 00270 // VALUE: Each row pair gives the 6x3 correct basis set values at an evaluation point: (P,F,D) layout 00271 double basisValues[] = { 00272 // 4 vertices 00273 1.0,0.,0., 0.,0.,0., 0.,-1.0,0., 0.,0.,1.0, 00274 0.,0.,0., 0.,0.,0., 00275 00276 1.0,1.0,1.0, 0.,1.,0., 0.,0.,0., 0.,0.,0., 00277 0.,0.,1., 0.,0.,0., 00278 00279 0.,0.,0., -1.,0.,0., -1.0,-1.0,-1.0, 00280 0.,0.,0., 0.,0.,0., 0.,0.,1., 00281 00282 0.,0.,0., 0.,0.,0., 0.,0.,0., 1.0,1.0,1.0, 00283 -1.,0.,0., 0.,-1.,0., 00284 00285 // 6 edge centers 00286 1.0,0.5,0.5, 0.,0.5,0., 0.,-0.5,0., 00287 0.,0.,0.5, 0.,0.,0.5, 0.,0.,0., 00288 00289 0.5,0.5,0.5, -0.5,0.5,0., 00290 -0.5,-0.5,-0.5, 0.,0.,0., 0.,0.,0.5, 0.,0.,0.5, 00291 00292 0.5,0.,0., -0.5,0.,0., -0.5,-1.0,-0.5, 00293 0.,0.,0.5, 0.,0.,0., 0.,0.,0.5, 00294 00295 0.5,0.,0., 0.,0.,0., 0.,-0.5,0., 0.5,0.5,1.0, 00296 -0.5,0.,0., 0.,-0.5,0., 00297 00298 0.5,0.5,0.5, 0.,0.5,0., 0.,0.,0., 0.5,0.5,0.5, 00299 -0.5,0.,0.5, 0.,-0.5,0., 00300 00301 0.,0.,0., -0.5,0.,0., -0.5,-0.5,-0.5, 0.5,0.5,0.5, 00302 -0.5,0.,0., 0.,-0.5,0.5 00303 }; 00304 00305 // CURL: each row pair gives the 3x12 correct values of the curls of the 12 basis functions: (P,F,D) layout 00306 double basisCurls[] = { 00307 // 4 vertices 00308 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00309 0.,-2.0,0., 2.0,0.,0., 00310 00311 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00312 0.,-2.0,0., 2.0,0.,0., 00313 00314 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00315 0.,-2.0,0., 2.0,0.,0., 00316 00317 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00318 0.,-2.0,0., 2.0,0.,0., 00319 00320 // 6 edge centers 00321 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00322 0.,-2.0,0., 2.0,0.,0., 00323 00324 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00325 0.,-2.0,0., 2.0,0.,0., 00326 00327 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00328 0.,-2.0,0., 2.0,0.,0., 00329 00330 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00331 0.,-2.0,0., 2.0,0.,0., 00332 00333 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00334 0.,-2.0,0., 2.0,0.,0., 00335 00336 0.,-2.0,2.0, 0.,0.,2.0, -2.0,0.,2.0, -2.0,2.0,0., 00337 0.,-2.0,0., 2.0,0.,0., 00338 }; 00339 00340 try{ 00341 00342 // Dimensions for the output arrays: 00343 int numFields = tetBasis.getCardinality(); 00344 int numPoints = tetNodes.dimension(0); 00345 int spaceDim = tetBasis.getBaseCellTopology().getDimension(); 00346 00347 // Generic array for values and curls that will be properly sized before each call 00348 FieldContainer<double> vals; 00349 00350 // Check VALUE of basis functions: resize vals to rank-3 container: 00351 vals.resize(numFields, numPoints, spaceDim); 00352 tetBasis.getValues(vals, tetNodes, OPERATOR_VALUE); 00353 for (int i = 0; i < numFields; i++) { 00354 for (int j = 0; j < numPoints; j++) { 00355 for (int k = 0; k < spaceDim; k++) { 00356 00357 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k 00358 int l = k + i * spaceDim + j * spaceDim * numFields; 00359 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) { 00360 errorFlag++; 00361 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00362 00363 // Output the multi-index of the value where the error is: 00364 *outStream << " At multi-index { "; 00365 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00366 *outStream << "} computed value: " << vals(i,j,k) 00367 << " but reference value: " << basisValues[l] << "\n"; 00368 } 00369 } 00370 } 00371 } 00372 00373 // Check CURL of basis function: resize vals to rank-3 container 00374 vals.resize(numFields, numPoints, spaceDim); 00375 tetBasis.getValues(vals, tetNodes, OPERATOR_CURL); 00376 for (int i = 0; i < numFields; i++) { 00377 for (int j = 0; j < numPoints; j++) { 00378 for (int k = 0; k < spaceDim; k++) { 00379 00380 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k 00381 int l = k + i * spaceDim + j * spaceDim * numFields; 00382 if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) { 00383 errorFlag++; 00384 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00385 00386 // Output the multi-index of the value where the error is: 00387 *outStream << " At multi-index { "; 00388 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00389 *outStream << "} computed curl component: " << vals(i,j,k) 00390 << " but reference curl component: " << basisCurls[l] << "\n"; 00391 } 00392 } 00393 } 00394 } 00395 00396 } 00397 00398 // Catch unexpected errors 00399 catch (std::logic_error err) { 00400 *outStream << err.what() << "\n\n"; 00401 errorFlag = -1000; 00402 }; 00403 00404 *outStream \ 00405 << "\n" 00406 << "===============================================================================\n"\ 00407 << "| TEST 4: correctness of DoF locations |\n"\ 00408 << "===============================================================================\n"; 00409 00410 try{ 00411 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis = 00412 Teuchos::rcp(new Basis_HCURL_TET_I1_FEM<double, FieldContainer<double> >); 00413 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface = 00414 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis); 00415 00416 int spaceDim = 3; 00417 FieldContainer<double> cvals; 00418 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality(),spaceDim); // last dimension is spatial dim 00419 00420 // Check exceptions. 00421 #ifdef HAVE_INTREPID_DEBUG 00422 cvals.resize(1,2,3); 00423 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00424 cvals.resize(3,2); 00425 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00426 cvals.resize(4,2); 00427 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00428 #endif 00429 cvals.resize(6,spaceDim); 00430 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--; 00431 // Check if number of thrown exceptions matches the one we expect 00432 if (throwCounter != nException) { 00433 errorFlag++; 00434 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00435 } 00436 00437 // Check mathematical correctness 00438 FieldContainer<double> tangents(basis->getCardinality(),spaceDim); // tangents at each point basis point 00439 tangents(0,0) = 1.0; tangents(0,1) = 0.0; tangents(0,2) = 0.0; 00440 tangents(1,0) = -1.0; tangents(1,1) = 1.0; tangents(1,2) = 0.0; 00441 tangents(2,0) = 0.0; tangents(2,1) = -1.0; tangents(2,2) = 0.0; 00442 tangents(3,0) = 0.0; tangents(3,1) = 0.0; tangents(3,2) = 1.0; 00443 tangents(4,0) = -1.0; tangents(4,1) = 0.0; tangents(4,2) = 1.0; 00444 tangents(5,0) = 0.0; tangents(5,1) = -1.0; tangents(5,2) = 1.0; 00445 00446 basis->getValues(bvals, cvals, OPERATOR_VALUE); 00447 char buffer[120]; 00448 for (int i=0; i<bvals.dimension(0); i++) { 00449 for (int j=0; j<bvals.dimension(1); j++) { 00450 00451 double tangent = 0.0; 00452 for(int d=0;d<spaceDim;d++) 00453 tangent += bvals(i,j,d)*tangents(j,d); 00454 00455 if ((i != j) && (std::abs(tangent - 0.0) > INTREPID_TOL)) { 00456 errorFlag++; 00457 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); 00458 *outStream << buffer; 00459 } 00460 else if ((i == j) && (std::abs(tangent - 1.0) > INTREPID_TOL)) { 00461 errorFlag++; 00462 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); 00463 *outStream << buffer; 00464 } 00465 } 00466 } 00467 00468 } 00469 catch (std::logic_error err){ 00470 *outStream << err.what() << "\n\n"; 00471 errorFlag = -1000; 00472 }; 00473 00474 if (errorFlag != 0) 00475 std::cout << "End Result: TEST FAILED\n"; 00476 else 00477 std::cout << "End Result: TEST PASSED\n"; 00478 00479 // reset format state of std::cout 00480 std::cout.copyfmt(oldFormatState); 00481 00482 return errorFlag; 00483 }
1.7.6.1