|
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_TRI_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_TRI_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_TRI_I1_FEM<double, FieldContainer<double> > triBasis; 00110 int errorFlag = 0; 00111 00112 // Define throw number for exception testing 00113 int nException = 0; 00114 int throwCounter = 0; 00115 00116 // Define array containing the 3 vertices of the reference TRI and its 3 edge midpoints. 00117 FieldContainer<double> triNodes(7, 2); 00118 triNodes(0,0) = 0.0; triNodes(0,1) = 0.0; 00119 triNodes(1,0) = 1.0; triNodes(1,1) = 0.0; 00120 triNodes(2,0) = 0.0; triNodes(2,1) = 1.0; 00121 // edge midpoints 00122 triNodes(3,0) = 0.5; triNodes(3,1) = 0.0; 00123 triNodes(4,0) = 0.5; triNodes(4,1) = 0.5; 00124 triNodes(5,0) = 0.0; triNodes(5,1) = 0.5; 00125 // Inside Triangle 00126 triNodes(6,0) = 0.25; triNodes(6,1) = 0.25; 00127 00128 // Generic array for the output values; needs to be properly resized depending on the operator type 00129 FieldContainer<double> vals; 00130 00131 00132 try{ 00133 // exception #1: GRAD cannot be applied to HCURL functions 00134 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00135 vals.resize(triBasis.getCardinality(), triNodes.dimension(0), 4 ); 00136 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_GRAD), throwCounter, nException ); 00137 00138 // exception #2: DIV cannot be applied to HCURL functions 00139 // resize vals to rank-2 container with dimensions (num. points, num. basis functions) 00140 vals.resize(triBasis.getCardinality(), triNodes.dimension(0) ); 00141 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_DIV), throwCounter, nException ); 00142 00143 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00144 // getDofTag() to access invalid array elements thereby causing bounds check exception 00145 // exception #3 00146 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00147 // exception #4 00148 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00149 // exception #5 00150 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(0,4,1), throwCounter, nException ); 00151 // exception #6 00152 INTREPID_TEST_COMMAND( triBasis.getDofTag(12), throwCounter, nException ); 00153 // exception #7 00154 INTREPID_TEST_COMMAND( triBasis.getDofTag(-1), throwCounter, nException ); 00155 00156 #ifdef HAVE_INTREPID_DEBUG 00157 // Exceptions 8-15 test exception handling with incorrectly dimensioned input/output arrays 00158 // exception #8: input points array must be of rank-2 00159 FieldContainer<double> badPoints1(4, 5, 3); 00160 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00161 00162 // exception #9 dimension 1 in the input point array must equal space dimension of the cell 00163 FieldContainer<double> badPoints2(4, triBasis.getBaseCellTopology().getDimension() + 1); 00164 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00165 00166 // exception #10 output values must be of rank-3 for OPERATOR_VALUE in 2D 00167 FieldContainer<double> badVals1(4, 3); 00168 INTREPID_TEST_COMMAND( triBasis.getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00169 00170 FieldContainer<double> badCurls1(4,3,2); 00171 // exception #11 output values must be of rank-2 for OPERATOR_CURL 00172 INTREPID_TEST_COMMAND( triBasis.getValues(badCurls1, triNodes, OPERATOR_CURL), throwCounter, nException ); 00173 00174 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions) 00175 FieldContainer<double> badVals2(triBasis.getCardinality() + 1, triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension()); 00176 INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00177 00178 // exception #13 incorrect 1st dimension of output array (must equal number of points) 00179 FieldContainer<double> badVals3(triBasis.getCardinality(), triNodes.dimension(0) + 1, triBasis.getBaseCellTopology().getDimension() ); 00180 INTREPID_TEST_COMMAND( triBasis.getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00181 00182 // exception #14: incorrect 2nd dimension of output array for VALUE (must equal the space dimension) 00183 FieldContainer<double> badVals4(triBasis.getCardinality(), triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension() - 1); 00184 INTREPID_TEST_COMMAND( triBasis.getValues(badVals4, triNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00185 00186 // exception #15: D2 cannot be applied to HCURL functions 00187 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00188 vals.resize(triBasis.getCardinality(), 00189 triNodes.dimension(0), 00190 Intrepid::getDkCardinality(OPERATOR_D2, triBasis.getBaseCellTopology().getDimension())); 00191 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_D2), 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 = triBasis.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 = triBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00221 00222 std::vector<int> myTag = triBasis.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 < triBasis.getCardinality(); bfOrd++) { 00244 std::vector<int> myTag = triBasis.getDofTag(bfOrd); 00245 int myBfOrd = triBasis.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: correct values in (P,F,D) layout 00275 double basisValues[] = { 00276 1.000, 0, 0, 0, 0, -1.000, 1.000, 1.000, 0, 1.000, 0, 0, 0, 0, \ 00277 -1.000, 0, -1.000, -1.000, 1.000, 0.5000, 0, 0.5000, 0, -0.5000, \ 00278 0.5000, 0.5000, -0.5000, 0.5000, -0.5000, -0.5000, 0.5000, 0, \ 00279 -0.5000, 0, -0.5000, -1.000, 0.7500, 0.2500, -0.2500, 0.2500, \ 00280 -0.2500, -0.7500}; 00281 00282 // CURL: correct values in (P,F) layout 00283 double basisCurls[] = { 00284 2.0, 2.0, 2.0, 00285 2.0, 2.0, 2.0, 00286 2.0, 2.0, 2.0, 00287 2.0, 2.0, 2.0, 00288 2.0, 2.0, 2.0, 00289 2.0, 2.0, 2.0, 00290 2.0, 2.0, 2.0 00291 }; 00292 00293 try{ 00294 00295 // Dimensions for the output arrays: 00296 int numFields = triBasis.getCardinality(); 00297 int numPoints = triNodes.dimension(0); 00298 int spaceDim = triBasis.getBaseCellTopology().getDimension(); 00299 00300 // Generic array for values and curls that will be properly sized before each call 00301 FieldContainer<double> vals; 00302 00303 // Check VALUE of basis functions: resize vals to rank-3 container: 00304 vals.resize(numFields, numPoints, spaceDim); 00305 triBasis.getValues(vals, triNodes, OPERATOR_VALUE); 00306 for (int i = 0; i < numFields; i++) { 00307 for (int j = 0; j < numPoints; j++) { 00308 for (int k = 0; k < spaceDim; k++) { 00309 00310 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k 00311 int l = k + i * spaceDim + j * spaceDim * numFields; 00312 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) { 00313 errorFlag++; 00314 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00315 00316 // Output the multi-index of the value where the error is: 00317 *outStream << " At multi-index { "; 00318 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00319 *outStream << "} computed value: " << vals(i,j,k) 00320 << " but reference value: " << basisValues[l] << "\n"; 00321 } 00322 } 00323 } 00324 } 00325 00326 // Check CURL of basis function: resize vals to rank-2 container 00327 vals.resize(numFields, numPoints); 00328 triBasis.getValues(vals, triNodes, OPERATOR_CURL); 00329 for (int i = 0; i < numFields; i++) { 00330 for (int j = 0; j < numPoints; j++) { 00331 int l = i + j * numFields; 00332 if (std::abs(vals(i,j) - basisCurls[l]) > INTREPID_TOL) { 00333 errorFlag++; 00334 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00335 00336 // Output the multi-index of the value where the error is: 00337 *outStream << " At multi-index { "; 00338 *outStream << i << " ";*outStream << j << " "; 00339 *outStream << "} computed curl component: " << vals(i,j) 00340 << " but reference curl component: " << basisCurls[l] << "\n"; 00341 } 00342 } 00343 } 00344 } 00345 00346 // Catch unexpected errors 00347 catch (std::logic_error err) { 00348 *outStream << err.what() << "\n\n"; 00349 errorFlag = -1000; 00350 }; 00351 00352 *outStream \ 00353 << "\n" 00354 << "===============================================================================\n"\ 00355 << "| TEST 4: correctness of DoF locations |\n"\ 00356 << "===============================================================================\n"; 00357 00358 try{ 00359 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis = 00360 Teuchos::rcp(new Basis_HCURL_TRI_I1_FEM<double, FieldContainer<double> >); 00361 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface = 00362 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis); 00363 00364 int spaceDim = 2; 00365 FieldContainer<double> cvals; 00366 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality(),2); // last dimension is spatial dim 00367 00368 // Check exceptions. 00369 #ifdef HAVE_INTREPID_DEBUG 00370 cvals.resize(1,2,3); 00371 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00372 cvals.resize(4,2); 00373 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00374 cvals.resize(4,3); 00375 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00376 #endif 00377 cvals.resize(3,spaceDim); 00378 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--; 00379 // Check if number of thrown exceptions matches the one we expect 00380 if (throwCounter != nException) { 00381 errorFlag++; 00382 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00383 } 00384 00385 // Check mathematical correctness 00386 FieldContainer<double> tangents(basis->getCardinality(),spaceDim); // tangents at each point basis point 00387 00388 tangents(0,0) = 1.0; tangents(0,1) = 0.0; 00389 tangents(1,0) = -1.0; tangents(1,1) = 1.0; 00390 tangents(2,0) = 0.0; tangents(2,1) = -1.0; 00391 00392 basis->getValues(bvals, cvals, OPERATOR_VALUE); 00393 char buffer[120]; 00394 for (int i=0; i<bvals.dimension(0); i++) { 00395 for (int j=0; j<bvals.dimension(1); j++) { 00396 00397 double tangent = 0.0; 00398 for(int d=0;d<spaceDim;d++) 00399 tangent += bvals(i,j,d)*tangents(j,d); 00400 00401 if ((i != j) && (std::abs(tangent - 0.0) > INTREPID_TOL)) { 00402 errorFlag++; 00403 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), tangent, 0.0); 00404 *outStream << buffer; 00405 } 00406 else if ((i == j) && (std::abs(tangent - 1.0) > INTREPID_TOL)) { 00407 errorFlag++; 00408 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), tangent, 1.0); 00409 *outStream << buffer; 00410 } 00411 } 00412 } 00413 00414 } 00415 catch (std::logic_error err){ 00416 *outStream << err.what() << "\n\n"; 00417 errorFlag = -1000; 00418 }; 00419 00420 00421 if (errorFlag != 0) 00422 std::cout << "End Result: TEST FAILED\n"; 00423 else 00424 std::cout << "End Result: TEST PASSED\n"; 00425 00426 // reset format state of std::cout 00427 std::cout.copyfmt(oldFormatState); 00428 00429 return errorFlag; 00430 }
1.7.6.1