|
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_HDIV_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_HDIV_TRI_I1_FEM) |\n" \ 00093 << "| |\n" \ 00094 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \ 00095 << "| 3) Exception tests on input arguments and invalid operators |\n" \ 00096 << "| 2) Basis values for VALUE and DIV operators |\n" \ 00097 << "| |\n" \ 00098 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00099 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00100 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00101 << "| |\n" \ 00102 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00103 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00104 << "| |\n" \ 00105 << "===============================================================================\n"\ 00106 << "| TEST 1: Basis creation, exception testing |\n"\ 00107 << "===============================================================================\n"; 00108 00109 // Define basis and error flag 00110 Basis_HDIV_TRI_I1_FEM<double, FieldContainer<double> > triBasis; 00111 int errorFlag = 0; 00112 00113 // Initialize throw counter for exception testing 00114 int nException = 0; 00115 int throwCounter = 0; 00116 00117 // Define array containing the 3 vertices of the reference TRI and its 3 edge midpoints. 00118 FieldContainer<double> triNodes(6, 2); 00119 triNodes(0,0) = 0.0; triNodes(0,1) = 0.0; 00120 triNodes(1,0) = 1.0; triNodes(1,1) = 0.0; 00121 triNodes(2,0) = 0.0; triNodes(2,1) = 1.0; 00122 // edge midpoints 00123 triNodes(3,0) = 0.5; triNodes(3,1) = 0.0; 00124 triNodes(4,0) = 0.5; triNodes(4,1) = 0.5; 00125 triNodes(5,0) = 0.0; triNodes(5,1) = 0.5; 00126 00127 00128 // Generic array for the output values; needs to be properly resized depending on the operator type 00129 FieldContainer<double> vals; 00130 00131 try{ 00132 // exception #1: GRAD cannot be applied to HDIV functions 00133 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00134 vals.resize(triBasis.getCardinality(), triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension() ); 00135 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_GRAD), throwCounter, nException ); 00136 00137 // exception #2: CURL cannot be applied to HDIV functions 00138 INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_CURL), throwCounter, nException ); 00139 00140 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00141 // getDofTag() to access invalid array elements thereby causing bounds check exception 00142 // exception #3 00143 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00144 // exception #4 00145 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00146 // exception #5 00147 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(0,2,1), throwCounter, nException ); 00148 // exception #6 00149 INTREPID_TEST_COMMAND( triBasis.getDofTag(3), throwCounter, nException ); 00150 // exception #7 00151 INTREPID_TEST_COMMAND( triBasis.getDofTag(-1), throwCounter, nException ); 00152 // exception #8 00153 INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(2,0,0), throwCounter, nException ); 00154 00155 #ifdef HAVE_INTREPID_DEBUG 00156 // Exceptions 9-16 test exception handling with incorrectly dimensioned input/output arrays 00157 // exception #9: input points array must be of rank-2 00158 FieldContainer<double> badPoints1(4, 5, 3); 00159 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00160 00161 // exception #10 dimension 1 in the input point array must equal space dimension of the cell 00162 FieldContainer<double> badPoints2(4, triBasis.getBaseCellTopology().getDimension() + 1); 00163 INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00164 00165 // exception #11 output values must be of rank-3 for OPERATOR_VALUE 00166 FieldContainer<double> badVals1(4, 3); 00167 INTREPID_TEST_COMMAND( triBasis.getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00168 00169 // exception #12 output values must be of rank-2 for OPERATOR_DIV 00170 FieldContainer<double> badVals2(4, 3, 1); 00171 INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00172 00173 // exception #13 incorrect 0th dimension of output array for OPERATOR_VALUE (must equal number of basis functions) 00174 FieldContainer<double> badVals3(triBasis.getCardinality() + 1, triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension()); 00175 INTREPID_TEST_COMMAND( triBasis.getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00176 00177 // exception #14 incorrect 0th dimension of output array for OPERATOR_DIV (must equal number of basis functions) 00178 FieldContainer<double> badVals4(triBasis.getCardinality() + 1, triNodes.dimension(0)); 00179 INTREPID_TEST_COMMAND( triBasis.getValues(badVals4, triNodes, OPERATOR_DIV), throwCounter, nException ); 00180 00181 // exception #15 incorrect 1st dimension of output array (must equal number of points) 00182 FieldContainer<double> badVals5(triBasis.getCardinality(), triNodes.dimension(0) + 1, triBasis.getBaseCellTopology().getDimension()); 00183 INTREPID_TEST_COMMAND( triBasis.getValues(badVals5, triNodes, OPERATOR_VALUE), throwCounter, nException ); 00184 00185 // exception #16 incorrect 1st dimension of output array (must equal number of points) 00186 FieldContainer<double> badVals6(triBasis.getCardinality(), triNodes.dimension(0) + 1); 00187 INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_DIV), throwCounter, nException ); 00188 00189 // exception #17: incorrect 2nd dimension of output array (must equal the space dimension) 00190 FieldContainer<double> badVals7(triBasis.getCardinality(), triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension() + 1); 00191 INTREPID_TEST_COMMAND( triBasis.getValues(badVals7, triNodes, OPERATOR_VALUE), 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 if (throwCounter != nException) { 00204 errorFlag++; 00205 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00206 } 00207 00208 *outStream \ 00209 << "\n" 00210 << "===============================================================================\n"\ 00211 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00212 << "===============================================================================\n"; 00213 00214 try{ 00215 std::vector<std::vector<int> > allTags = triBasis.getAllDofTags(); 00216 00217 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00218 for (unsigned i = 0; i < allTags.size(); i++) { 00219 int bfOrd = triBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00220 00221 std::vector<int> myTag = triBasis.getDofTag(bfOrd); 00222 if( !( (myTag[0] == allTags[i][0]) && 00223 (myTag[1] == allTags[i][1]) && 00224 (myTag[2] == allTags[i][2]) && 00225 (myTag[3] == allTags[i][3]) ) ) { 00226 errorFlag++; 00227 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00228 *outStream << " getDofOrdinal( {" 00229 << allTags[i][0] << ", " 00230 << allTags[i][1] << ", " 00231 << allTags[i][2] << ", " 00232 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00233 *outStream << " getDofTag(" << bfOrd << ") = { " 00234 << myTag[0] << ", " 00235 << myTag[1] << ", " 00236 << myTag[2] << ", " 00237 << myTag[3] << "}\n"; 00238 } 00239 } 00240 00241 // Now do the same but loop over basis functions 00242 for( int bfOrd = 0; bfOrd < triBasis.getCardinality(); bfOrd++) { 00243 std::vector<int> myTag = triBasis.getDofTag(bfOrd); 00244 int myBfOrd = triBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00245 if( bfOrd != myBfOrd) { 00246 errorFlag++; 00247 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00248 *outStream << " getDofTag(" << bfOrd << ") = { " 00249 << myTag[0] << ", " 00250 << myTag[1] << ", " 00251 << myTag[2] << ", " 00252 << myTag[3] << "} but getDofOrdinal({" 00253 << myTag[0] << ", " 00254 << myTag[1] << ", " 00255 << myTag[2] << ", " 00256 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00257 } 00258 } 00259 } 00260 catch (std::logic_error err){ 00261 *outStream << err.what() << "\n\n"; 00262 errorFlag = -1000; 00263 }; 00264 00265 *outStream \ 00266 << "\n" 00267 << "===============================================================================\n"\ 00268 << "| TEST 3: correctness of basis function values |\n"\ 00269 << "===============================================================================\n"; 00270 00271 outStream -> precision(20); 00272 00273 // VALUE: Correct basis values in (F,P,D) format: each group of two rows gives basis function 00274 // values at vertices followed by midpoints. This is the same array format as the output from getValues. 00275 double basisValues[] = { 00276 // basis function 0 at 3 vertices followed by 3 midpoints 00277 0.0,-1.0, 1.0,-1.0, 0.0, 0.0, 00278 0.5,-1.0, 0.5,-0.5, 0.0,-0.5, 00279 // basis function 1 at 3 vertices followed by 3 midpoints 00280 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 00281 0.5, 0.0, 0.5, 0.5, 0.0, 0.5, 00282 // basis function 2 at 3 vertices followed by 3 midpoints 00283 -1.0, 0.0, 0.0, 0.0, -1.0, 1.0, 00284 -0.5, 0.0, -0.5, 0.5, -1.0, 0.5 00285 }; 00286 00287 // DIV: each row gives the 3 correct values of the divergence of the 3 basis functions 00288 double basisDivs[] = { 00289 // 3 vertices 00290 2.0, 2.0, 2.0, 00291 2.0, 2.0, 2.0, 00292 2.0, 2.0, 2.0, 00293 // 3 edge centers 00294 2.0, 2.0, 2.0, 00295 2.0, 2.0, 2.0, 00296 2.0, 2.0, 2.0, 00297 }; 00298 00299 try{ 00300 00301 // Dimensions for the output arrays: 00302 int numFields = triBasis.getCardinality(); 00303 int numPoints = triNodes.dimension(0); 00304 int spaceDim = triBasis.getBaseCellTopology().getDimension(); 00305 00306 // Generic array for values and divs that will be properly sized before each call 00307 FieldContainer<double> vals; 00308 00309 // Check VALUE of basis functions: resize vals to rank-3 container: 00310 vals.resize(numFields, numPoints, spaceDim); 00311 triBasis.getValues(vals, triNodes, OPERATOR_VALUE); 00312 for (int i = 0; i < numFields; i++) { 00313 for (int j = 0; j < numPoints; j++) { 00314 for (int k = 0; k < spaceDim; k++) { 00315 // basisValues are in (F,P,D) format and the multiindex is (i,j,k), here's the offset: 00316 int l = k + j * spaceDim + i * spaceDim * numPoints; 00317 00318 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) { 00319 errorFlag++; 00320 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00321 00322 // Output the multi-index of the value where the error is: 00323 *outStream << " address = "<< l <<"\n"; 00324 *outStream << " At multi-index { "; 00325 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00326 *outStream << "} computed value: " << vals(i,j,k) 00327 << " but reference value: " << basisValues[l] << "\n"; 00328 } 00329 } 00330 } 00331 } 00332 00333 00334 // Check DIV of basis function: resize vals to rank-2 container 00335 vals.resize(numFields, numPoints); 00336 triBasis.getValues(vals, triNodes, OPERATOR_DIV); 00337 for (int i = 0; i < numFields; i++) { 00338 for (int j = 0; j < numPoints; j++) { 00339 int l = i + j * numFields; 00340 if (std::abs(vals(i,j) - basisDivs[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 (Field,Point,Dim) multi-index { "; 00346 *outStream << i << " ";*outStream << j << " "; 00347 *outStream << "} computed divergence component: " << vals(i,j) 00348 << " but reference divergence component: " << basisDivs[l] << "\n"; 00349 } 00350 } 00351 } 00352 00353 } 00354 00355 // Catch unexpected errors 00356 catch (std::logic_error err) { 00357 *outStream << err.what() << "\n\n"; 00358 errorFlag = -1000; 00359 }; 00360 00361 if (errorFlag != 0) 00362 std::cout << "End Result: TEST FAILED\n"; 00363 else 00364 std::cout << "End Result: TEST PASSED\n"; 00365 00366 // reset format state of std::cout 00367 std::cout.copyfmt(oldFormatState); 00368 00369 return errorFlag; 00370 }
1.7.6.1