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