|
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_QUAD_In_FEM.hpp" 00050 #include "Intrepid_PointTools.hpp" 00051 #include "Teuchos_oblackholestream.hpp" 00052 #include "Teuchos_RCP.hpp" 00053 #include "Teuchos_GlobalMPISession.hpp" 00054 00055 using namespace std; 00056 using namespace Intrepid; 00057 00058 #define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \ 00059 { \ 00060 ++nException; \ 00061 try { \ 00062 S ; \ 00063 } \ 00064 catch (std::logic_error err) { \ 00065 ++throwCounter; \ 00066 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \ 00067 *outStream << err.what() << '\n'; \ 00068 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \ 00069 }; \ 00070 } 00071 00072 int main(int argc, char *argv[]) { 00073 00074 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00075 00076 // This little trick lets us print to std::cout only if 00077 // a (dummy) command-line argument is provided. 00078 int iprint = argc - 1; 00079 Teuchos::RCP<std::ostream> outStream; 00080 Teuchos::oblackholestream bhs; // outputs nothing 00081 if (iprint > 0) 00082 outStream = Teuchos::rcp(&std::cout, false); 00083 else 00084 outStream = Teuchos::rcp(&bhs, false); 00085 00086 // Save the format state of the original std::cout. 00087 Teuchos::oblackholestream oldFormatState; 00088 oldFormatState.copyfmt(std::cout); 00089 00090 *outStream \ 00091 << "===============================================================================\n" \ 00092 << "| |\n" \ 00093 << "| Unit Test (Basis_HCURL_QUAD_In_FEM) |\n" \ 00094 << "| |\n" \ 00095 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \ 00096 << "| 2) Basis values for VALUE and CURL 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 const int deg = 1; 00111 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >()); 00112 00113 Basis_HCURL_QUAD_In_FEM<double, FieldContainer<double> > quadBasis(deg,POINTTYPE_EQUISPACED); 00114 int errorFlag = 0; 00115 00116 // Initialize throw counter for exception testing 00117 int nException = 0; 00118 int throwCounter = 0; 00119 00120 // Array with the 4 vertices of the reference Quadrilateral, its center and 4 more points 00121 FieldContainer<double> quadNodes(9, 2); 00122 quadNodes(0,0) = -1.0; quadNodes(0,1) = -1.0; 00123 quadNodes(1,0) = 0.0; quadNodes(1,1) = -1.0; 00124 quadNodes(2,0) = 1.0; quadNodes(2,1) = -1.0; 00125 quadNodes(3,0) = -1.0; quadNodes(3,1) = 0.0; 00126 quadNodes(4,0) = 0.0; quadNodes(4,1) = 0.0; 00127 quadNodes(5,0) = 1.0; quadNodes(5,1) = 0.0; 00128 quadNodes(6,0) = -1.0; quadNodes(6,1) = 1.0; 00129 quadNodes(7,0) = 0.0; quadNodes(7,1) = 1.0; 00130 quadNodes(8,0) = 1.0; quadNodes(8,1) = 1.0; 00131 00132 00133 // Generic array for the output values; needs to be properly resized depending on the operator type 00134 FieldContainer<double> vals; 00135 00136 try{ 00137 // exception #1: GRAD cannot be applied to HCURL functions 00138 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00139 vals.resize(quadBasis.getCardinality(), quadNodes.dimension(0), 4 ); 00140 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_GRAD), throwCounter, nException ); 00141 00142 // exception #2: DIV cannot be applied to HCURL functions 00143 // resize vals to rank-2 container with dimensions (num. points, num. basis functions) 00144 vals.resize(quadBasis.getCardinality(), quadNodes.dimension(0) ); 00145 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_DIV), throwCounter, nException ); 00146 00147 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00148 // getDofTag() to access invalid array elements thereby causing bounds check exception 00149 // exception #3 00150 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00151 // exception #4 00152 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00153 // exception #5 00154 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(0,4,1), throwCounter, nException ); 00155 // exception #6 00156 INTREPID_TEST_COMMAND( quadBasis.getDofTag(12), throwCounter, nException ); 00157 // exception #7 00158 INTREPID_TEST_COMMAND( quadBasis.getDofTag(-1), throwCounter, nException ); 00159 00160 #ifdef HAVE_INTREPID_DEBUG 00161 // Exceptions 8-15 test exception handling with incorrectly dimensioned input/output arrays 00162 // exception #8: input points array must be of rank-2 00163 FieldContainer<double> badPoints1(4, 5, 3); 00164 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00165 00166 // exception #9 dimension 1 in the input point array must equal space dimension of the cell 00167 FieldContainer<double> badPoints2(4, quadBasis.getBaseCellTopology().getDimension() + 1); 00168 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00169 00170 // exception #10 output values must be of rank-3 for OPERATOR_VALUE in 2D 00171 FieldContainer<double> badVals1(4, 3); 00172 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals1, quadNodes, OPERATOR_VALUE), throwCounter, nException ); 00173 00174 FieldContainer<double> badCurls1(4,3,2); 00175 // exception #11 output values must be of rank-2 for OPERATOR_CURL 00176 INTREPID_TEST_COMMAND( quadBasis.getValues(badCurls1, quadNodes, OPERATOR_CURL), throwCounter, nException ); 00177 00178 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions) 00179 FieldContainer<double> badVals2(quadBasis.getCardinality() + 1, quadNodes.dimension(0), quadBasis.getBaseCellTopology().getDimension()); 00180 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00181 00182 // exception #13 incorrect 1st dimension of output array (must equal number of points) 00183 FieldContainer<double> badVals3(quadBasis.getCardinality(), quadNodes.dimension(0) + 1, quadBasis.getBaseCellTopology().getDimension() ); 00184 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals3, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00185 00186 // exception #14: incorrect 2nd dimension of output array for VALUE (must equal the space dimension) 00187 FieldContainer<double> badVals4(quadBasis.getCardinality(), quadNodes.dimension(0), quadBasis.getBaseCellTopology().getDimension() - 1); 00188 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals4, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ; 00189 00190 // exception #15: D2 cannot be applied to HCURL functions 00191 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary) 00192 vals.resize(quadBasis.getCardinality(), 00193 quadNodes.dimension(0), 00194 Intrepid::getDkCardinality(OPERATOR_D2, quadBasis.getBaseCellTopology().getDimension())); 00195 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_D2), throwCounter, nException ); 00196 #endif 00197 00198 } 00199 catch (std::logic_error err) { 00200 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00201 *outStream << err.what() << '\n'; 00202 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00203 errorFlag = -1000; 00204 }; 00205 00206 // Check if number of thrown exceptions matches the one we expect 00207 // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match. 00208 if (throwCounter != nException) { 00209 errorFlag++; 00210 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00211 } 00212 //#endif 00213 00214 *outStream \ 00215 << "\n" 00216 << "===============================================================================\n"\ 00217 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00218 << "===============================================================================\n"; 00219 00220 try{ 00221 std::vector<std::vector<int> > allTags = quadBasis.getAllDofTags(); 00222 00223 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00224 for (unsigned i = 0; i < allTags.size(); i++) { 00225 int bfOrd = quadBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00226 00227 std::vector<int> myTag = quadBasis.getDofTag(bfOrd); 00228 if( !( (myTag[0] == allTags[i][0]) && 00229 (myTag[1] == allTags[i][1]) && 00230 (myTag[2] == allTags[i][2]) && 00231 (myTag[3] == allTags[i][3]) ) ) { 00232 errorFlag++; 00233 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00234 *outStream << " getDofOrdinal( {" 00235 << allTags[i][0] << ", " 00236 << allTags[i][1] << ", " 00237 << allTags[i][2] << ", " 00238 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00239 *outStream << " getDofTag(" << bfOrd << ") = { " 00240 << myTag[0] << ", " 00241 << myTag[1] << ", " 00242 << myTag[2] << ", " 00243 << myTag[3] << "}\n"; 00244 } 00245 } 00246 00247 // Now do the same but loop over basis functions 00248 for( int bfOrd = 0; bfOrd < quadBasis.getCardinality(); bfOrd++) { 00249 std::vector<int> myTag = quadBasis.getDofTag(bfOrd); 00250 int myBfOrd = quadBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00251 if( bfOrd != myBfOrd) { 00252 errorFlag++; 00253 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00254 *outStream << " getDofTag(" << bfOrd << ") = { " 00255 << myTag[0] << ", " 00256 << myTag[1] << ", " 00257 << myTag[2] << ", " 00258 << myTag[3] << "} but getDofOrdinal({" 00259 << myTag[0] << ", " 00260 << myTag[1] << ", " 00261 << myTag[2] << ", " 00262 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00263 } 00264 } 00265 } 00266 catch (std::logic_error err){ 00267 *outStream << err.what() << "\n\n"; 00268 errorFlag = -1000; 00269 }; 00270 00271 *outStream \ 00272 << "\n" 00273 << "===============================================================================\n"\ 00274 << "| TEST 3: correctness of basis function values |\n"\ 00275 << "===============================================================================\n"; 00276 00277 outStream -> precision(20); 00278 00279 // VALUE: Each row pair gives the 4x2 correct basis set values at an evaluation point: (P,F,D) layout 00280 double basisValues[] = { 00281 // first bf, first row of points 00282 1,0,1,0,1,0, 00283 // first bf, second row of points 00284 0.5,0,0.5,0,0.5,0, 00285 // first bf, third row of points 00286 0,0,0,0,0,0, 00287 // second bf, first row of points 00288 0,0,0,0,0,0, 00289 // second bf, second row of points 00290 0.5,0,0.5,0,0.5,0, 00291 // second bf, third row of points 00292 1,0,1,0,1,0, 00293 // third bf, first row of points 00294 0,1,0,0.5,0,0, 00295 // third bf, second row of points 00296 0,1,0,0.5,0,0, 00297 // third bf, third row of points 00298 0,1,0,0.5,0,0, 00299 // fourth bf, first row of points 00300 0,0,0,0.5,0,1, 00301 // fourth bf, second row of points 00302 0,0,0,0.5,0,1, 00303 // fourth bf, third row of points 00304 0,0,0,0.5,0,1, 00305 }; 00306 00307 // CURL: correct values in (F,P) format 00308 double basisCurls[] = { 00309 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 00310 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, 00311 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, 00312 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5 00313 }; 00314 00315 00316 try{ 00317 // Dimensions for the output arrays: 00318 int numFields = quadBasis.getCardinality(); 00319 int numPoints = quadNodes.dimension(0); 00320 int spaceDim = quadBasis.getBaseCellTopology().getDimension(); 00321 00322 // Generic array for values and curls that will be properly sized before each call 00323 FieldContainer<double> vals; 00324 00325 // Check VALUE of basis functions: resize vals to rank-3 container: 00326 vals.resize(numFields, numPoints, spaceDim); 00327 quadBasis.getValues(vals, quadNodes, OPERATOR_VALUE); 00328 for (int i = 0; i < numFields; i++) { 00329 for (int j = 0; j < numPoints; j++) { 00330 for (int k = 0; k < spaceDim; k++) { 00331 00332 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k 00333 int l = k + i * spaceDim * numPoints + j * spaceDim; 00334 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) { 00335 errorFlag++; 00336 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00337 00338 // Output the multi-index of the value where the error is: 00339 *outStream << " At multi-index { "; 00340 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00341 *outStream << "} computed value: " << vals(i,j,k) 00342 << " but reference value: " << basisValues[l] << "\n"; 00343 } 00344 } 00345 } 00346 } 00347 00348 // Check CURL of basis function: resize vals to rank-2 container 00349 vals.resize(numFields, numPoints); 00350 quadBasis.getValues(vals, quadNodes, OPERATOR_CURL); 00351 for (int i = 0; i < numFields; i++) { 00352 for (int j = 0; j < numPoints; j++) { 00353 int l = j + i * numPoints; 00354 if (std::abs(vals(i,j) - basisCurls[l]) > INTREPID_TOL) { 00355 errorFlag++; 00356 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00357 00358 // Output the multi-index of the value where the error is: 00359 *outStream << " At multi-index { "; 00360 *outStream << i << " ";*outStream << j << " "; 00361 *outStream << "} computed curl component: " << vals(i,j) 00362 << " but reference curl component: " << basisCurls[l] << "\n"; 00363 } 00364 } 00365 } 00366 00367 } 00368 00369 // Catch unexpected errors 00370 catch (std::logic_error err) { 00371 *outStream << err.what() << "\n\n"; 00372 errorFlag = -1000; 00373 }; 00374 00375 if (errorFlag != 0) 00376 std::cout << "End Result: TEST FAILED\n"; 00377 else 00378 std::cout << "End Result: TEST PASSED\n"; 00379 00380 // reset format state of std::cout 00381 std::cout.copyfmt(oldFormatState); 00382 00383 return errorFlag; 00384 }
1.7.6.1