|
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 // This library is free software; you can redistribute it and/or modify 00011 // it under the terms of the GNU Lesser General Public License as 00012 // published by the Free Software Foundation; either version 2.1 of the 00013 // License, or (at your option) any later version. 00014 // 00015 // This library is distributed in the hope that it will be useful, but 00016 // WITHOUT ANY WARRANTY; without even the implied warranty of 00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00018 // Lesser General Public License for more details. 00019 // 00020 // You should have received a copy of the GNU Lesser General Public 00021 // License along with this library; if not, write to the Free Software 00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 00023 // USA 00024 // Questions? Contact Pavel Bochev (pbboche@sandia.gov), 00025 // Denis Ridzal (dridzal@sandia.gov), 00026 // Kara Peterson (kjpeter@sandia.gov). 00027 // 00028 // ************************************************************************ 00029 // @HEADER 00030 00035 #include "Intrepid_FieldContainer.hpp" 00036 #include "Intrepid_HGRAD_QUAD_C1_FEM.hpp" 00037 #include "Teuchos_oblackholestream.hpp" 00038 #include "Teuchos_RCP.hpp" 00039 #include "Teuchos_GlobalMPISession.hpp" 00040 00041 using namespace std; 00042 using namespace Intrepid; 00043 00044 #define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \ 00045 { \ 00046 ++nException; \ 00047 try { \ 00048 S ; \ 00049 } \ 00050 catch (std::logic_error err) { \ 00051 ++throwCounter; \ 00052 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \ 00053 *outStream << err.what() << '\n'; \ 00054 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \ 00055 }; \ 00056 } 00057 00058 int main(int argc, char *argv[]) { 00059 00060 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00061 00062 // This little trick lets us print to std::cout only if 00063 // a (dummy) command-line argument is provided. 00064 int iprint = argc - 1; 00065 Teuchos::RCP<std::ostream> outStream; 00066 Teuchos::oblackholestream bhs; // outputs nothing 00067 if (iprint > 0) 00068 outStream = Teuchos::rcp(&std::cout, false); 00069 else 00070 outStream = Teuchos::rcp(&bhs, false); 00071 00072 // Save the format state of the original std::cout. 00073 Teuchos::oblackholestream oldFormatState; 00074 oldFormatState.copyfmt(std::cout); 00075 00076 *outStream \ 00077 << "===============================================================================\n" \ 00078 << "| |\n" \ 00079 << "| Unit Test (Basis_HGRAD_QUAD_C1_FEM) |\n" \ 00080 << "| |\n" \ 00081 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \ 00082 << "| 2) Basis values for VALUE, GRAD, CURL, and Dk operators |\n" \ 00083 << "| |\n" \ 00084 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00085 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00086 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00087 << "| |\n" \ 00088 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00089 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00090 << "| |\n" \ 00091 << "===============================================================================\n"\ 00092 << "| TEST 1: Basis creation, exception testing |\n"\ 00093 << "===============================================================================\n"; 00094 00095 // Define basis and error flag 00096 Basis_HGRAD_QUAD_C1_FEM<double, FieldContainer<double> > quadBasis; 00097 int errorFlag = 0; 00098 00099 // Initialize throw counter for exception testing 00100 int nException = 0; 00101 int throwCounter = 0; 00102 00103 // Define array containing the 4 vertices of the reference QUAD and its center. 00104 FieldContainer<double> quadNodes(5, 2); 00105 quadNodes(0,0) = -1.0; quadNodes(0,1) = -1.0; 00106 quadNodes(1,0) = 1.0; quadNodes(1,1) = -1.0; 00107 quadNodes(2,0) = 1.0; quadNodes(2,1) = 1.0; 00108 quadNodes(3,0) = -1.0; quadNodes(3,1) = 1.0; 00109 quadNodes(4,0) = 0.0; quadNodes(4,1) = 0.0; 00110 00111 // Generic array for the output values; needs to be properly resized depending on the operator type 00112 FieldContainer<double> vals; 00113 00114 try{ 00115 // exception #1: DIV cannot be applied to scalar functions 00116 // resize vals to rank-2 container with dimensions (num. points, num. basis functions) 00117 vals.resize(quadBasis.getCardinality(), quadNodes.dimension(0)); 00118 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_DIV), throwCounter, nException ); 00119 00120 // Exceptions 2-6: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00121 // getDofTag() to access invalid array elements thereby causing bounds check exception 00122 // exception #2 00123 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00124 // exception #3 00125 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00126 // exception #4 00127 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(0,4,0), throwCounter, nException ); 00128 // exception #5 00129 INTREPID_TEST_COMMAND( quadBasis.getDofTag(5), throwCounter, nException ); 00130 // exception #6 00131 INTREPID_TEST_COMMAND( quadBasis.getDofTag(-1), throwCounter, nException ); 00132 00133 #ifdef HAVE_INTREPID_DEBUG 00134 // Exceptions 7- test exception handling with incorrectly dimensioned input/output arrays 00135 // exception #7: input points array must be of rank-2 00136 FieldContainer<double> badPoints1(4, 5, 3); 00137 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00138 00139 // exception #8 dimension 1 in the input point array must equal space dimension of the cell 00140 FieldContainer<double> badPoints2(4, 3); 00141 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00142 00143 // exception #9 output values must be of rank-2 for OPERATOR_VALUE 00144 FieldContainer<double> badVals1(4, 3, 1); 00145 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals1, quadNodes, OPERATOR_VALUE), throwCounter, nException ); 00146 00147 // exception #10 output values must be of rank-3 for OPERATOR_GRAD 00148 FieldContainer<double> badVals2(4, 3); 00149 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_GRAD), throwCounter, nException ); 00150 00151 // exception #11 output values must be of rank-3 for OPERATOR_CURL 00152 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_CURL), throwCounter, nException ); 00153 00154 // exception #12 output values must be of rank-3 for OPERATOR_D2 00155 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_D2), throwCounter, nException ); 00156 00157 // exception #13 incorrect 0th dimension of output array (must equal number of basis functions) 00158 FieldContainer<double> badVals3(quadBasis.getCardinality() + 1, quadNodes.dimension(0)); 00159 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals3, quadNodes, OPERATOR_VALUE), throwCounter, nException ); 00160 00161 // exception #14 incorrect 1st dimension of output array (must equal number of points) 00162 FieldContainer<double> badVals4(quadBasis.getCardinality(), quadNodes.dimension(0) + 1); 00163 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals4, quadNodes, OPERATOR_VALUE), throwCounter, nException ); 00164 00165 // exception #15: incorrect 2nd dimension of output array (must equal the space dimension) 00166 FieldContainer<double> badVals5(quadBasis.getCardinality(), quadNodes.dimension(0), 4); 00167 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals5, quadNodes, OPERATOR_GRAD), throwCounter, nException ); 00168 00169 // exception #16: incorrect 2nd dimension of output array (must equal D2 cardinality in 2D) 00170 FieldContainer<double> badVals6(quadBasis.getCardinality(), quadNodes.dimension(0), 40); 00171 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals6, quadNodes, OPERATOR_D2), throwCounter, nException ); 00172 00173 // exception #17: incorrect 2nd dimension of output array (must equal D3 cardinality in 2D) 00174 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals6, quadNodes, OPERATOR_D3), throwCounter, nException ); 00175 #endif 00176 00177 } 00178 catch (std::logic_error err) { 00179 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00180 *outStream << err.what() << '\n'; 00181 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00182 errorFlag = -1000; 00183 }; 00184 00185 // Check if number of thrown exceptions matches the one we expect 00186 if (throwCounter != nException) { 00187 errorFlag++; 00188 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00189 } 00190 00191 *outStream \ 00192 << "\n" 00193 << "===============================================================================\n"\ 00194 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00195 << "===============================================================================\n"; 00196 00197 try{ 00198 std::vector<std::vector<int> > allTags = quadBasis.getAllDofTags(); 00199 00200 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00201 for (unsigned i = 0; i < allTags.size(); i++) { 00202 int bfOrd = quadBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00203 00204 std::vector<int> myTag = quadBasis.getDofTag(bfOrd); 00205 if( !( (myTag[0] == allTags[i][0]) && 00206 (myTag[1] == allTags[i][1]) && 00207 (myTag[2] == allTags[i][2]) && 00208 (myTag[3] == allTags[i][3]) ) ) { 00209 errorFlag++; 00210 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00211 *outStream << " getDofOrdinal( {" 00212 << allTags[i][0] << ", " 00213 << allTags[i][1] << ", " 00214 << allTags[i][2] << ", " 00215 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00216 *outStream << " getDofTag(" << bfOrd << ") = { " 00217 << myTag[0] << ", " 00218 << myTag[1] << ", " 00219 << myTag[2] << ", " 00220 << myTag[3] << "}\n"; 00221 } 00222 } 00223 00224 // Now do the same but loop over basis functions 00225 for( int bfOrd = 0; bfOrd < quadBasis.getCardinality(); bfOrd++) { 00226 std::vector<int> myTag = quadBasis.getDofTag(bfOrd); 00227 int myBfOrd = quadBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00228 if( bfOrd != myBfOrd) { 00229 errorFlag++; 00230 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00231 *outStream << " getDofTag(" << bfOrd << ") = { " 00232 << myTag[0] << ", " 00233 << myTag[1] << ", " 00234 << myTag[2] << ", " 00235 << myTag[3] << "} but getDofOrdinal({" 00236 << myTag[0] << ", " 00237 << myTag[1] << ", " 00238 << myTag[2] << ", " 00239 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00240 } 00241 } 00242 } 00243 catch (std::logic_error err){ 00244 *outStream << err.what() << "\n\n"; 00245 errorFlag = -1000; 00246 }; 00247 00248 *outStream \ 00249 << "\n" 00250 << "===============================================================================\n"\ 00251 << "| TEST 3: correctness of basis function values |\n"\ 00252 << "===============================================================================\n"; 00253 00254 outStream -> precision(20); 00255 00256 // VALUE: Each row gives the 4 correct basis set values at an evaluation point 00257 double basisValues[] = { 00258 1.0, 0.0, 0.0, 0.0, 00259 0.0, 1.0, 0.0, 0.0, 00260 0.0, 0.0, 1.0, 0.0, 00261 0.0, 0.0, 0.0, 1.0, 00262 0.25,0.25,0.25,0.25 00263 }; 00264 00265 // GRAD and D1: each row gives the 8 correct values of the gradients of the 4 basis functions 00266 double basisGrads[] = { 00267 -0.5, -0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.5, 00268 -0.5, 0.0, 0.5, -0.5, 0.0, 0.5, 0.0, 0.0, 00269 0.0, 0.0, 0.0, -0.5, 0.5, 0.5, -0.5, 0.0, 00270 0.0, -0.5, 0.0, 0.0, 0.5, 0.0, -0.5, 0.5, 00271 -0.25,-0.25, 0.25,-0.25, 0.25, 0.25, -0.25, 0.25 00272 }; 00273 00274 // CURL: each row gives the 8 correct values of the curls of the 4 basis functions 00275 double basisCurls[] = { 00276 -0.5, 0.5, 0.0, -0.5, 0.0, 0.0, 0.5, 0.0, 00277 0.0, 0.5, -0.5, -0.5, 0.5, 0.0, 0.0, 0.0, 00278 0.0, 0.0, -0.5, 0.0, 0.5, -0.5, 0.0, 0.5, 00279 -0.5, 0.0, 0.0, 0.0, 0.0, -0.5, 0.5, 0.5, 00280 -0.25, 0.25, -0.25,-0.25, 0.25,-0.25, 0.25, 0.25 00281 }; 00282 00283 //D2: each row gives the 12 correct values of all 2nd derivatives of the 4 basis functions 00284 double basisD2[] = { 00285 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 00286 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 00287 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 00288 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 00289 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 0.0, 0.25, 0.0, 0.0,-0.25, 0.0, 00290 }; 00291 00292 try{ 00293 00294 // Dimensions for the output arrays: 00295 int numFields = quadBasis.getCardinality(); 00296 int numPoints = quadNodes.dimension(0); 00297 int spaceDim = quadBasis.getBaseCellTopology().getDimension(); 00298 int D2Cardin = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim); 00299 00300 // Generic array for values, grads, curls, etc. that will be properly sized before each call 00301 FieldContainer<double> vals; 00302 00303 // Check VALUE of basis functions: resize vals to rank-2 container: 00304 vals.resize(numFields, numPoints); 00305 quadBasis.getValues(vals, quadNodes, OPERATOR_VALUE); 00306 for (int i = 0; i < numFields; i++) { 00307 for (int j = 0; j < numPoints; j++) { 00308 int l = i + j * numFields; 00309 if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) { 00310 errorFlag++; 00311 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00312 00313 // Output the multi-index of the value where the error is: 00314 *outStream << " At multi-index { "; 00315 *outStream << i << " ";*outStream << j << " "; 00316 *outStream << "} computed value: " << vals(i,j) 00317 << " but reference value: " << basisValues[l] << "\n"; 00318 } 00319 } 00320 } 00321 00322 // Check GRAD of basis function: resize vals to rank-3 container 00323 vals.resize(numFields, numPoints, spaceDim); 00324 quadBasis.getValues(vals, quadNodes, OPERATOR_GRAD); 00325 for (int i = 0; i < numFields; i++) { 00326 for (int j = 0; j < numPoints; j++) { 00327 for (int k = 0; k < spaceDim; k++) { 00328 int l = k + i * spaceDim + j * spaceDim * numFields; 00329 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00330 errorFlag++; 00331 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00332 00333 // Output the multi-index of the value where the error is: 00334 *outStream << " At multi-index { "; 00335 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00336 *outStream << "} computed grad component: " << vals(i,j,k) 00337 << " but reference grad component: " << basisGrads[l] << "\n"; 00338 } 00339 } 00340 } 00341 } 00342 00343 00344 // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD) 00345 quadBasis.getValues(vals, quadNodes, OPERATOR_D1); 00346 for (int i = 0; i < numFields; i++) { 00347 for (int j = 0; j < numPoints; j++) { 00348 for (int k = 0; k < spaceDim; k++) { 00349 int l = k + i * spaceDim + j * spaceDim * numFields; 00350 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00351 errorFlag++; 00352 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00353 00354 // Output the multi-index of the value where the error is: 00355 *outStream << " At multi-index { "; 00356 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00357 *outStream << "} computed D1 component: " << vals(i,j,k) 00358 << " but reference D1 component: " << basisGrads[l] << "\n"; 00359 } 00360 } 00361 } 00362 } 00363 00364 00365 // Check CURL of basis function: resize vals just for illustration! 00366 vals.resize(numFields, numPoints, spaceDim); 00367 quadBasis.getValues(vals, quadNodes, OPERATOR_CURL); 00368 for (int i = 0; i < numFields; i++) { 00369 for (int j = 0; j < numPoints; j++) { 00370 for (int k = 0; k < spaceDim; k++) { 00371 int l = k + i * spaceDim + j * spaceDim * numFields; 00372 if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) { 00373 errorFlag++; 00374 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00375 00376 // Output the multi-index of the value where the error is: 00377 *outStream << " At multi-index { "; 00378 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00379 *outStream << "} computed curl component: " << vals(i,j,k) 00380 << " but reference curl component: " << basisCurls[l] << "\n"; 00381 } 00382 } 00383 } 00384 } 00385 00386 // Check D2 of basis function 00387 vals.resize(numFields, numPoints, D2Cardin); 00388 quadBasis.getValues(vals, quadNodes, OPERATOR_D2); 00389 for (int i = 0; i < numFields; i++) { 00390 for (int j = 0; j < numPoints; j++) { 00391 for (int k = 0; k < D2Cardin; k++) { 00392 int l = k + i * D2Cardin + j * D2Cardin * numFields; 00393 if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) { 00394 errorFlag++; 00395 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00396 00397 // Output the multi-index of the value where the error is: 00398 *outStream << " At multi-index { "; 00399 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00400 *outStream << "} computed D2 component: " << vals(i,j,k) 00401 << " but reference D2 component: " << basisD2[l] << "\n"; 00402 } 00403 } 00404 } 00405 } 00406 00407 00408 // Check all higher derivatives - must be zero. 00409 for(EOperator op = OPERATOR_D3; op < OPERATOR_MAX; op++) { 00410 00411 // The last dimension is the number of kth derivatives and needs to be resized for every Dk 00412 int DkCardin = Intrepid::getDkCardinality(op, spaceDim); 00413 vals.resize(numFields, numPoints, DkCardin); 00414 00415 quadBasis.getValues(vals, quadNodes, op); 00416 for (int i = 0; i < vals.size(); i++) { 00417 if (std::abs(vals[i]) > INTREPID_TOL) { 00418 errorFlag++; 00419 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00420 00421 // Get the multi-index of the value where the error is and the operator order 00422 std::vector<int> myIndex; 00423 vals.getMultiIndex(myIndex,i); 00424 int ord = Intrepid::getOperatorOrder(op); 00425 *outStream << " At multi-index { "; 00426 for(int j = 0; j < vals.rank(); j++) { 00427 *outStream << myIndex[j] << " "; 00428 } 00429 *outStream << "} computed D"<< ord <<" component: " << vals[i] 00430 << " but reference D" << ord << " component: 0 \n"; 00431 } 00432 } 00433 } 00434 } 00435 // Catch unexpected errors 00436 catch (std::logic_error err) { 00437 *outStream << err.what() << "\n\n"; 00438 errorFlag = -1000; 00439 }; 00440 00441 00442 *outStream \ 00443 << "\n" 00444 << "===============================================================================\n"\ 00445 << "| TEST 4: correctness of DoF locations |\n"\ 00446 << "===============================================================================\n"; 00447 00448 try{ 00449 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis = 00450 Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM<double, FieldContainer<double> >); 00451 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface = 00452 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis); 00453 00454 FieldContainer<double> cvals; 00455 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality()); 00456 00457 // Check exceptions. 00458 #ifdef HAVE_INTREPID_DEBUG 00459 cvals.resize(1,2,3); 00460 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00461 cvals.resize(3,2); 00462 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00463 cvals.resize(4,3); 00464 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00465 #endif 00466 cvals.resize(4,2); 00467 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--; 00468 // Check if number of thrown exceptions matches the one we expect 00469 if (throwCounter != nException) { 00470 errorFlag++; 00471 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00472 } 00473 00474 // Check mathematical correctness. 00475 basis->getValues(bvals, cvals, OPERATOR_VALUE); 00476 char buffer[120]; 00477 for (int i=0; i<bvals.dimension(0); i++) { 00478 for (int j=0; j<bvals.dimension(1); j++) { 00479 if ((i != j) && (std::abs(bvals(i,j) - 0.0) > INTREPID_TOL)) { 00480 errorFlag++; 00481 sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), bvals(i,j), 0.0); 00482 *outStream << buffer; 00483 } 00484 else if ((i == j) && (std::abs(bvals(i,j) - 1.0) > INTREPID_TOL)) { 00485 errorFlag++; 00486 sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), bvals(i,j), 1.0); 00487 *outStream << buffer; 00488 } 00489 } 00490 } 00491 00492 } 00493 catch (std::logic_error err){ 00494 *outStream << err.what() << "\n\n"; 00495 errorFlag = -1000; 00496 }; 00497 00498 if (errorFlag != 0) 00499 std::cout << "End Result: TEST FAILED\n"; 00500 else 00501 std::cout << "End Result: TEST PASSED\n"; 00502 00503 // reset format state of std::cout 00504 std::cout.copyfmt(oldFormatState); 00505 00506 return errorFlag; 00507 }
1.7.4