|
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 00036 #include "Intrepid_FieldContainer.hpp" 00037 #include "Intrepid_HGRAD_TRI_C1_FEM.hpp" 00038 #include "Teuchos_oblackholestream.hpp" 00039 #include "Teuchos_RCP.hpp" 00040 #include "Teuchos_GlobalMPISession.hpp" 00041 00042 using namespace std; 00043 using namespace Intrepid; 00044 00045 #define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \ 00046 { \ 00047 ++nException; \ 00048 try { \ 00049 S ; \ 00050 } \ 00051 catch (std::logic_error err) { \ 00052 ++throwCounter; \ 00053 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \ 00054 *outStream << err.what() << '\n'; \ 00055 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \ 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_TRI_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 (dridzal@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 00096 // Define basis and error flag 00097 Basis_HGRAD_TRI_C1_FEM<double, FieldContainer<double> > triBasis; 00098 int errorFlag = 0; 00099 00100 // Initialize throw counter for exception testing 00101 int nException = 0; 00102 int throwCounter = 0; 00103 00104 // Define array containing the 3 vertices of the reference Triangle, its center and another point 00105 FieldContainer<double> triNodes(5, 2); 00106 triNodes(0,0) = 0.0; triNodes(0,1) = 0.0; 00107 triNodes(1,0) = 1.0; triNodes(1,1) = 0.0; 00108 triNodes(2,0) = 0.0; triNodes(2,1) = 1.0; 00109 triNodes(3,0) = 0.5; triNodes(3,1) = 0.5; 00110 triNodes(4,0) = 0.0; triNodes(4,1) = 0.75; 00111 00112 // Generic array for the output values; needs to be properly resized depending on the operator type 00113 FieldContainer<double> vals; 00114 00115 try{ 00116 // exception #1: DIV cannot be applied to scalar functions 00117 // resize vals to rank-2 container with dimensions (num. points, num. basis functions) 00118 vals.resize(triBasis.getCardinality(), triNodes.dimension(0) ); 00119 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_DIV), throwCounter, nException ); 00120 00121 // Exceptions 2-6: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00122 // getDofTag() to access invalid array elements thereby causing bounds check exception 00123 // exception #2 00124 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(2,0,0), throwCounter, nException ); 00125 // exception #3 00126 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00127 // exception #4 00128 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(0,4,0), throwCounter, nException ); 00129 // exception #5 00130 INTREPID_TEST_COMMAND( triBasis.getDofTag(5), throwCounter, nException ); 00131 // exception #6 00132 INTREPID_TEST_COMMAND( triBasis.getDofTag(-1), throwCounter, nException ); 00133 00134 #ifdef HAVE_INTREPID_DEBUG 00135 // Exceptions 7-17 test exception handling with incorrectly dimensioned input/output arrays 00136 // exception #7: input points array must be of rank-2 00137 FieldContainer<double> badPoints1(4, 5, 3); 00138 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00139 00140 // exception #8 dimension 1 in the input point array must equal space dimension of the cell 00141 FieldContainer<double> badPoints2(4, 3); 00142 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00143 00144 // exception #9 output values must be of rank-2 for OPERATOR_VALUE 00145 FieldContainer<double> badVals1(4, 3, 1); 00146 INTREPID_TEST_COMMAND( triBasis.getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00147 00148 // exception #10 output values must be of rank-3 for OPERATOR_GRAD 00149 FieldContainer<double> badVals2(4, 3); 00150 INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_GRAD), throwCounter, nException ); 00151 00152 // exception #11 output values must be of rank-3 for OPERATOR_CURL 00153 INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_CURL), throwCounter, nException ); 00154 00155 // exception #12 output values must be of rank-3 for OPERATOR_D2 00156 INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_D2), throwCounter, nException ); 00157 00158 // exception #13 incorrect 1st dimension of output array (must equal number of basis functions) 00159 FieldContainer<double> badVals3(triBasis.getCardinality() + 1, triNodes.dimension(0)); 00160 INTREPID_TEST_COMMAND( triBasis.getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00161 00162 // exception #14 incorrect 0th dimension of output array (must equal number of points) 00163 FieldContainer<double> badVals4(triBasis.getCardinality(), triNodes.dimension(0) + 1); 00164 INTREPID_TEST_COMMAND( triBasis.getValues(badVals4, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00165 00166 // exception #15: incorrect 2nd dimension of output array (must equal the space dimension) 00167 FieldContainer<double> badVals5(triBasis.getCardinality(), triNodes.dimension(0), 4); 00168 INTREPID_TEST_COMMAND( triBasis.getValues(badVals5, triNodes, OPERATOR_GRAD), throwCounter, nException ); 00169 00170 // exception #16: incorrect 2nd dimension of output array (must equal D2 cardinality in 2D) 00171 FieldContainer<double> badVals6(triBasis.getCardinality(), triNodes.dimension(0), 40); 00172 INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_D2), throwCounter, nException ); 00173 00174 // exception #17: incorrect 2nd dimension of output array (must equal D3 cardinality in 2D) 00175 INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_D3), throwCounter, nException ); 00176 #endif 00177 00178 } 00179 catch (std::logic_error err) { 00180 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00181 *outStream << err.what() << '\n'; 00182 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00183 errorFlag = -1000; 00184 }; 00185 00186 // Check if number of thrown exceptions matches the one we expect 00187 if (throwCounter != nException) { 00188 errorFlag++; 00189 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00190 } 00191 00192 *outStream \ 00193 << "\n" 00194 << "===============================================================================\n"\ 00195 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00196 << "===============================================================================\n"; 00197 00198 try{ 00199 std::vector<std::vector<int> > allTags = triBasis.getAllDofTags(); 00200 00201 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00202 for (unsigned i = 0; i < allTags.size(); i++) { 00203 int bfOrd = triBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00204 00205 std::vector<int> myTag = triBasis.getDofTag(bfOrd); 00206 if( !( (myTag[0] == allTags[i][0]) && 00207 (myTag[1] == allTags[i][1]) && 00208 (myTag[2] == allTags[i][2]) && 00209 (myTag[3] == allTags[i][3]) ) ) { 00210 errorFlag++; 00211 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00212 *outStream << " getDofOrdinal( {" 00213 << allTags[i][0] << ", " 00214 << allTags[i][1] << ", " 00215 << allTags[i][2] << ", " 00216 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00217 *outStream << " getDofTag(" << bfOrd << ") = { " 00218 << myTag[0] << ", " 00219 << myTag[1] << ", " 00220 << myTag[2] << ", " 00221 << myTag[3] << "}\n"; 00222 } 00223 } 00224 00225 // Now do the same but loop over basis functions 00226 for( int bfOrd = 0; bfOrd < triBasis.getCardinality(); bfOrd++) { 00227 std::vector<int> myTag = triBasis.getDofTag(bfOrd); 00228 int myBfOrd = triBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00229 if( bfOrd != myBfOrd) { 00230 errorFlag++; 00231 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00232 *outStream << " getDofTag(" << bfOrd << ") = { " 00233 << myTag[0] << ", " 00234 << myTag[1] << ", " 00235 << myTag[2] << ", " 00236 << myTag[3] << "} but getDofOrdinal({" 00237 << myTag[0] << ", " 00238 << myTag[1] << ", " 00239 << myTag[2] << ", " 00240 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00241 } 00242 } 00243 } 00244 catch (std::logic_error err){ 00245 *outStream << err.what() << "\n\n"; 00246 errorFlag = -1000; 00247 }; 00248 00249 *outStream \ 00250 << "\n" 00251 << "===============================================================================\n"\ 00252 << "| TEST 3: correctness of basis function values |\n"\ 00253 << "===============================================================================\n"; 00254 00255 outStream -> precision(20); 00256 00257 // VALUE: Each row gives the 3 correct basis set values at an evaluation point 00258 double basisValues[] = { 00259 1.0, 0.0, 0.0, 00260 0.0, 1.0, 0.0, 00261 0.0, 0.0, 1.0, 00262 0.0, 0.5, 0.5, 00263 0.25,0.0, 0.75 00264 }; 00265 00266 // GRAD and D1: each row gives the 6 correct values of the gradients of the 3 basis functions 00267 double basisGrads[] = { 00268 -1.0, -1.0, 1.0, 0.0, 0.0, 1.0, 00269 -1.0, -1.0, 1.0, 0.0, 0.0, 1.0, 00270 -1.0, -1.0, 1.0, 0.0, 0.0, 1.0, 00271 -1.0, -1.0, 1.0, 0.0, 0.0, 1.0, 00272 -1.0, -1.0, 1.0, 0.0, 0.0, 1.0, 00273 }; 00274 00275 // CURL: each row gives the 6 correct values of the curls of the 3 basis functions 00276 double basisCurls[] = { 00277 -1.0, 1.0, 0.0, -1.0, 1.0, 0.0, 00278 -1.0, 1.0, 0.0, -1.0, 1.0, 0.0, 00279 -1.0, 1.0, 0.0, -1.0, 1.0, 0.0, 00280 -1.0, 1.0, 0.0, -1.0, 1.0, 0.0, 00281 -1.0, 1.0, 0.0, -1.0, 1.0, 0.0 00282 }; 00283 00284 try{ 00285 00286 // Dimensions for the output arrays: 00287 int numFields = triBasis.getCardinality(); 00288 int numPoints = triNodes.dimension(0); 00289 int spaceDim = triBasis.getBaseCellTopology().getDimension(); 00290 00291 // Generic array for values, grads, curls, etc. that will be properly sized before each call 00292 FieldContainer<double> vals; 00293 00294 // Check VALUE of basis functions: resize vals to rank-2 container: 00295 vals.resize(numFields, numPoints); 00296 triBasis.getValues(vals, triNodes, OPERATOR_VALUE); 00297 for (int i = 0; i < numFields; i++) { 00298 for (int j = 0; j < numPoints; j++) { 00299 int l = i + j * numFields; 00300 if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) { 00301 errorFlag++; 00302 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00303 00304 // Output the multi-index of the value where the error is: 00305 *outStream << " At multi-index { "; 00306 *outStream << i << " ";*outStream << j << " "; 00307 *outStream << "} computed value: " << vals(i,j) 00308 << " but reference value: " << basisValues[l] << "\n"; 00309 } 00310 } 00311 } 00312 00313 // Check GRAD of basis function: resize vals to rank-3 container 00314 vals.resize(numFields, numPoints, spaceDim); 00315 triBasis.getValues(vals, triNodes, OPERATOR_GRAD); 00316 for (int i = 0; i < numFields; i++) { 00317 for (int j = 0; j < numPoints; j++) { 00318 for (int k = 0; k < spaceDim; k++) { 00319 int l = k + i * spaceDim + j * spaceDim * numFields; 00320 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00321 errorFlag++; 00322 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00323 00324 // Output the multi-index of the value where the error is: 00325 *outStream << " At multi-index { "; 00326 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00327 *outStream << "} computed grad component: " << vals(i,j,k) 00328 << " but reference grad component: " << basisGrads[l] << "\n"; 00329 } 00330 } 00331 } 00332 } 00333 00334 // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD) 00335 triBasis.getValues(vals, triNodes, OPERATOR_D1); 00336 for (int i = 0; i < numFields; i++) { 00337 for (int j = 0; j < numPoints; j++) { 00338 for (int k = 0; k < spaceDim; k++) { 00339 int l = k + i * spaceDim + j * spaceDim * numFields; 00340 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00341 errorFlag++; 00342 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00343 00344 // Output the multi-index of the value where the error is: 00345 *outStream << " At multi-index { "; 00346 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00347 *outStream << "} computed D1 component: " << vals(i,j,k) 00348 << " but reference D1 component: " << basisGrads[l] << "\n"; 00349 } 00350 } 00351 } 00352 } 00353 00354 00355 // Check CURL of basis function: resize vals just for illustration! 00356 vals.resize(numFields, numPoints, spaceDim); 00357 triBasis.getValues(vals, triNodes, OPERATOR_CURL); 00358 for (int i = 0; i < numFields; i++) { 00359 for (int j = 0; j < numPoints; j++) { 00360 for (int k = 0; k < spaceDim; k++) { 00361 int l = k + i * spaceDim + j * spaceDim * numFields; 00362 if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) { 00363 errorFlag++; 00364 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00365 00366 // Output the multi-index of the value where the error is: 00367 *outStream << " At multi-index { "; 00368 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00369 *outStream << "} computed curl component: " << vals(i,j,k) 00370 << " but reference curl component: " << basisCurls[l] << "\n"; 00371 } 00372 } 00373 } 00374 } 00375 00376 // Check all higher derivatives - must be zero. 00377 for(EOperator op = OPERATOR_D2; op < OPERATOR_MAX; op++) { 00378 00379 // The last dimension is the number of kth derivatives and needs to be resized for every Dk 00380 int DkCardin = Intrepid::getDkCardinality(op, spaceDim); 00381 vals.resize(numFields, numPoints, DkCardin); 00382 00383 triBasis.getValues(vals, triNodes, op); 00384 for (int i = 0; i < vals.size(); i++) { 00385 if (std::abs(vals[i]) > INTREPID_TOL) { 00386 errorFlag++; 00387 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00388 00389 // Get the multi-index of the value where the error is and the operator order 00390 std::vector<int> myIndex; 00391 vals.getMultiIndex(myIndex,i); 00392 int ord = Intrepid::getOperatorOrder(op); 00393 *outStream << " At multi-index { "; 00394 for(int j = 0; j < vals.rank(); j++) { 00395 *outStream << myIndex[j] << " "; 00396 } 00397 *outStream << "} computed D"<< ord <<" component: " << vals[i] 00398 << " but reference D" << ord << " component: 0 \n"; 00399 } 00400 } 00401 } 00402 } 00403 00404 // Catch unexpected errors 00405 catch (std::logic_error err) { 00406 *outStream << err.what() << "\n\n"; 00407 errorFlag = -1000; 00408 }; 00409 00410 if (errorFlag != 0) 00411 std::cout << "End Result: TEST FAILED\n"; 00412 else 00413 std::cout << "End Result: TEST PASSED\n"; 00414 00415 // reset format state of std::cout 00416 std::cout.copyfmt(oldFormatState); 00417 00418 return errorFlag; 00419 }
1.7.4