|
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_HCURL_TRI_I1_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_HCURL_TRI_I1_FEM) |\n" \ 00080 << "| |\n" \ 00081 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \ 00082 << "| 2) Basis values for VALUE and CURL 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_HCURL_TRI_I1_FEM<double, FieldContainer<double> > triBasis; 00097 int errorFlag = 0; 00098 00099 // Define throw number for exception testing 00100 int nException = 0; 00101 int throwCounter = 0; 00102 00103 // Define array containing the 3 vertices of the reference TRI and its 3 edge midpoints. 00104 FieldContainer<double> triNodes(7, 2); 00105 triNodes(0,0) = 0.0; triNodes(0,1) = 0.0; 00106 triNodes(1,0) = 1.0; triNodes(1,1) = 0.0; 00107 triNodes(2,0) = 0.0; triNodes(2,1) = 1.0; 00108 // edge midpoints 00109 triNodes(3,0) = 0.5; triNodes(3,1) = 0.0; 00110 triNodes(4,0) = 0.5; triNodes(4,1) = 0.5; 00111 triNodes(5,0) = 0.0; triNodes(5,1) = 0.5; 00112 // Inside Triangle 00113 triNodes(6,0) = 0.25; triNodes(6,1) = 0.25; 00114 00115 // Generic array for the output values; needs to be properly resized depending on the operator type 00116 FieldContainer<double> vals; 00117 00118 00119 try{ 00120 // exception #1: GRAD cannot be applied to HCURL functions 00121 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00122 vals.resize(triBasis.getCardinality(), triNodes.dimension(0), 4 ); 00123 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_GRAD), throwCounter, nException ); 00124 00125 // exception #2: DIV cannot be applied to HCURL functions 00126 // resize vals to rank-2 container with dimensions (num. points, num. basis functions) 00127 vals.resize(triBasis.getCardinality(), triNodes.dimension(0) ); 00128 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_DIV), throwCounter, nException ); 00129 00130 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00131 // getDofTag() to access invalid array elements thereby causing bounds check exception 00132 // exception #3 00133 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00134 // exception #4 00135 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00136 // exception #5 00137 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(0,4,1), throwCounter, nException ); 00138 // exception #6 00139 INTREPID_TEST_COMMAND( triBasis.getDofTag(12), throwCounter, nException ); 00140 // exception #7 00141 INTREPID_TEST_COMMAND( triBasis.getDofTag(-1), throwCounter, nException ); 00142 00143 #ifdef HAVE_INTREPID_DEBUG 00144 // Exceptions 8-15 test exception handling with incorrectly dimensioned input/output arrays 00145 // exception #8: input points array must be of rank-2 00146 FieldContainer<double> badPoints1(4, 5, 3); 00147 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00148 00149 // exception #9 dimension 1 in the input point array must equal space dimension of the cell 00150 FieldContainer<double> badPoints2(4, triBasis.getBaseCellTopology().getDimension() + 1); 00151 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00152 00153 // exception #10 output values must be of rank-3 for OPERATOR_VALUE in 2D 00154 FieldContainer<double> badVals1(4, 3); 00155 INTREPID_TEST_COMMAND( triBasis.getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00156 00157 FieldContainer<double> badCurls1(4,3,2); 00158 // exception #11 output values must be of rank-2 for OPERATOR_CURL 00159 INTREPID_TEST_COMMAND( triBasis.getValues(badCurls1, triNodes, OPERATOR_CURL), throwCounter, nException ); 00160 00161 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions) 00162 FieldContainer<double> badVals2(triBasis.getCardinality() + 1, triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension()); 00163 INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00164 00165 // exception #13 incorrect 1st dimension of output array (must equal number of points) 00166 FieldContainer<double> badVals3(triBasis.getCardinality(), triNodes.dimension(0) + 1, triBasis.getBaseCellTopology().getDimension() ); 00167 INTREPID_TEST_COMMAND( triBasis.getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00168 00169 // exception #14: incorrect 2nd dimension of output array for VALUE (must equal the space dimension) 00170 FieldContainer<double> badVals4(triBasis.getCardinality(), triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension() - 1); 00171 INTREPID_TEST_COMMAND( triBasis.getValues(badVals4, triNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00172 00173 // exception #15: D2 cannot be applied to HCURL functions 00174 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00175 vals.resize(triBasis.getCardinality(), 00176 triNodes.dimension(0), 00177 Intrepid::getDkCardinality(OPERATOR_D2, triBasis.getBaseCellTopology().getDimension())); 00178 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_D2), throwCounter, nException ); 00179 #endif 00180 00181 } 00182 catch (std::logic_error err) { 00183 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00184 *outStream << err.what() << '\n'; 00185 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00186 errorFlag = -1000; 00187 }; 00188 00189 // Check if number of thrown exceptions matches the one we expect 00190 // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match. 00191 if (throwCounter != nException) { 00192 errorFlag++; 00193 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00194 } 00195 00196 *outStream \ 00197 << "\n" 00198 << "===============================================================================\n"\ 00199 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00200 << "===============================================================================\n"; 00201 00202 try{ 00203 std::vector<std::vector<int> > allTags = triBasis.getAllDofTags(); 00204 00205 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00206 for (unsigned i = 0; i < allTags.size(); i++) { 00207 int bfOrd = triBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00208 00209 std::vector<int> myTag = triBasis.getDofTag(bfOrd); 00210 if( !( (myTag[0] == allTags[i][0]) && 00211 (myTag[1] == allTags[i][1]) && 00212 (myTag[2] == allTags[i][2]) && 00213 (myTag[3] == allTags[i][3]) ) ) { 00214 errorFlag++; 00215 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00216 *outStream << " getDofOrdinal( {" 00217 << allTags[i][0] << ", " 00218 << allTags[i][1] << ", " 00219 << allTags[i][2] << ", " 00220 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00221 *outStream << " getDofTag(" << bfOrd << ") = { " 00222 << myTag[0] << ", " 00223 << myTag[1] << ", " 00224 << myTag[2] << ", " 00225 << myTag[3] << "}\n"; 00226 } 00227 } 00228 00229 // Now do the same but loop over basis functions 00230 for( int bfOrd = 0; bfOrd < triBasis.getCardinality(); bfOrd++) { 00231 std::vector<int> myTag = triBasis.getDofTag(bfOrd); 00232 int myBfOrd = triBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00233 if( bfOrd != myBfOrd) { 00234 errorFlag++; 00235 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00236 *outStream << " getDofTag(" << bfOrd << ") = { " 00237 << myTag[0] << ", " 00238 << myTag[1] << ", " 00239 << myTag[2] << ", " 00240 << myTag[3] << "} but getDofOrdinal({" 00241 << myTag[0] << ", " 00242 << myTag[1] << ", " 00243 << myTag[2] << ", " 00244 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00245 } 00246 } 00247 } 00248 catch (std::logic_error err){ 00249 *outStream << err.what() << "\n\n"; 00250 errorFlag = -1000; 00251 }; 00252 00253 *outStream \ 00254 << "\n" 00255 << "===============================================================================\n"\ 00256 << "| TEST 3: correctness of basis function values |\n"\ 00257 << "===============================================================================\n"; 00258 00259 outStream -> precision(20); 00260 00261 // VALUE: correct values in (P,F,D) layout 00262 double basisValues[] = { 00263 1.000, 0, 0, 0, 0, -1.000, 1.000, 1.000, 0, 1.000, 0, 0, 0, 0, \ 00264 -1.000, 0, -1.000, -1.000, 1.000, 0.5000, 0, 0.5000, 0, -0.5000, \ 00265 0.5000, 0.5000, -0.5000, 0.5000, -0.5000, -0.5000, 0.5000, 0, \ 00266 -0.5000, 0, -0.5000, -1.000, 0.7500, 0.2500, -0.2500, 0.2500, \ 00267 -0.2500, -0.7500}; 00268 00269 // CURL: correct values in (P,F) layout 00270 double basisCurls[] = { 00271 2.0, 2.0, 2.0, 00272 2.0, 2.0, 2.0, 00273 2.0, 2.0, 2.0, 00274 2.0, 2.0, 2.0, 00275 2.0, 2.0, 2.0, 00276 2.0, 2.0, 2.0, 00277 2.0, 2.0, 2.0 00278 }; 00279 00280 try{ 00281 00282 // Dimensions for the output arrays: 00283 int numFields = triBasis.getCardinality(); 00284 int numPoints = triNodes.dimension(0); 00285 int spaceDim = triBasis.getBaseCellTopology().getDimension(); 00286 00287 // Generic array for values and curls that will be properly sized before each call 00288 FieldContainer<double> vals; 00289 00290 // Check VALUE of basis functions: resize vals to rank-3 container: 00291 vals.resize(numFields, numPoints, spaceDim); 00292 triBasis.getValues(vals, triNodes, OPERATOR_VALUE); 00293 for (int i = 0; i < numFields; i++) { 00294 for (int j = 0; j < numPoints; j++) { 00295 for (int k = 0; k < spaceDim; k++) { 00296 00297 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k 00298 int l = k + i * spaceDim + j * spaceDim * numFields; 00299 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) { 00300 errorFlag++; 00301 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00302 00303 // Output the multi-index of the value where the error is: 00304 *outStream << " At multi-index { "; 00305 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00306 *outStream << "} computed value: " << vals(i,j,k) 00307 << " but reference value: " << basisValues[l] << "\n"; 00308 } 00309 } 00310 } 00311 } 00312 00313 // Check CURL of basis function: resize vals to rank-2 container 00314 vals.resize(numFields, numPoints); 00315 triBasis.getValues(vals, triNodes, OPERATOR_CURL); 00316 for (int i = 0; i < numFields; i++) { 00317 for (int j = 0; j < numPoints; j++) { 00318 int l = i + j * numFields; 00319 if (std::abs(vals(i,j) - basisCurls[l]) > INTREPID_TOL) { 00320 errorFlag++; 00321 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00322 00323 // Output the multi-index of the value where the error is: 00324 *outStream << " At multi-index { "; 00325 *outStream << i << " ";*outStream << j << " "; 00326 *outStream << "} computed curl component: " << vals(i,j) 00327 << " but reference curl component: " << basisCurls[l] << "\n"; 00328 } 00329 } 00330 } 00331 } 00332 00333 00334 // Catch unexpected errors 00335 catch (std::logic_error err) { 00336 *outStream << err.what() << "\n\n"; 00337 errorFlag = -1000; 00338 }; 00339 00340 if (errorFlag != 0) 00341 std::cout << "End Result: TEST FAILED\n"; 00342 else 00343 std::cout << "End Result: TEST PASSED\n"; 00344 00345 // reset format state of std::cout 00346 std::cout.copyfmt(oldFormatState); 00347 00348 return errorFlag; 00349 }
1.7.4