|
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_HGRAD_QUAD_C2_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_HGRAD_QUAD_C2_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 (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_HGRAD_QUAD_C2_FEM<double, FieldContainer<double> > quadBasis; 00097 int errorFlag = 0; 00098 00099 // Initialize throw counter for exception testing 00100 int nException = 0; 00101 int throwCounter = 0; 00102 00103 // Array with the 4 vertices, 4 edge midpoints, center of the reference QUAD and a random point. 00104 FieldContainer<double> quadNodes(10, 2); 00105 quadNodes(0,0) = -1.0; quadNodes(0,1) = -1.0; 00106 quadNodes(1,0) = 1.0; quadNodes(1,1) = -1.0; 00107 quadNodes(2,0) = 1.0; quadNodes(2,1) = 1.0; 00108 quadNodes(3,0) = -1.0; quadNodes(3,1) = 1.0; 00109 // edge midpoints 00110 quadNodes(4,0) = 0.0; quadNodes(4,1) = -1.0; 00111 quadNodes(5,0) = 1.0; quadNodes(5,1) = 0.0; 00112 quadNodes(6,0) = 0.0; quadNodes(6,1) = 1.0; 00113 quadNodes(7,0) = -1.0; quadNodes(7,1) = 0.0; 00114 // center & random point 00115 quadNodes(8,0) = 0.0; quadNodes(8,1) = 0.0; 00116 quadNodes(9,0) =1./3.; quadNodes(9,1) =-3./5.; 00117 00118 00119 // Generic array for the output values; needs to be properly resized depending on the operator type 00120 FieldContainer<double> vals; 00121 00122 try{ 00123 // exception #1: DIV cannot be applied to scalar functions 00124 // resize vals to rank-2 container with dimensions (num. points, num. basis functions) 00125 vals.resize(quadBasis.getCardinality(), quadNodes.dimension(0)); 00126 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_DIV), throwCounter, nException ); 00127 00128 // Exceptions 2-6: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and 00129 // getDofTag() to access invalid array elements thereby causing bounds check exception 00130 // exception #2 00131 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(3,0,0), throwCounter, nException ); 00132 // exception #3 00133 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(1,1,1), throwCounter, nException ); 00134 // exception #4 00135 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(0,4,0), throwCounter, nException ); 00136 // exception #5 00137 INTREPID_TEST_COMMAND( quadBasis.getDofTag(10), throwCounter, nException ); 00138 // exception #6 00139 INTREPID_TEST_COMMAND( quadBasis.getDofTag(-1), throwCounter, nException ); 00140 00141 #ifdef HAVE_INTREPID_DEBUG 00142 // Exceptions 7- test exception handling with incorrectly dimensioned input/output arrays 00143 // exception #7: input points array must be of rank-2 00144 FieldContainer<double> badPoints1(4, 5, 3); 00145 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException ); 00146 00147 // exception #8 dimension 1 in the input point array must equal space dimension of the cell 00148 FieldContainer<double> badPoints2(4, quadBasis.getBaseCellTopology().getDimension() + 1); 00149 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException ); 00150 00151 // exception #9 output values must be of rank-2 for OPERATOR_VALUE 00152 FieldContainer<double> badVals1(4, 3, 1); 00153 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals1, quadNodes, OPERATOR_VALUE), throwCounter, nException ); 00154 00155 // exception #10 output values must be of rank-3 for OPERATOR_GRAD 00156 FieldContainer<double> badVals2(4, 3); 00157 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_GRAD), throwCounter, nException ); 00158 00159 // exception #11 output values must be of rank-3 for OPERATOR_CURL 00160 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_CURL), throwCounter, nException ); 00161 00162 // exception #12 output values must be of rank-3 for OPERATOR_D2 00163 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_D2), throwCounter, nException ); 00164 00165 // exception #13 incorrect 0th dimension of output array (must equal number of basis functions) 00166 FieldContainer<double> badVals3(quadBasis.getCardinality() + 1, quadNodes.dimension(0)); 00167 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals3, quadNodes, OPERATOR_VALUE), throwCounter, nException ); 00168 00169 // exception #14 incorrect 1st dimension of output array (must equal number of points in quadNodes) 00170 FieldContainer<double> badVals4(quadBasis.getCardinality(), quadNodes.dimension(0) + 1); 00171 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals4, quadNodes, OPERATOR_VALUE), throwCounter, nException ); 00172 00173 // exception #15: incorrect 2nd dimension of output array (must equal the space dimension) 00174 FieldContainer<double> badVals5(quadBasis.getCardinality(), quadNodes.dimension(0), quadBasis.getBaseCellTopology().getDimension() + 1); 00175 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals5, quadNodes, OPERATOR_GRAD), throwCounter, nException ); 00176 00177 // exception #16: incorrect 2nd dimension of output array (must equal D2 cardinality in 2D) 00178 FieldContainer<double> badVals6(quadBasis.getCardinality(), quadNodes.dimension(0), 40); 00179 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals6, quadNodes, OPERATOR_D2), throwCounter, nException ); 00180 00181 // exception #17: incorrect 2nd dimension of output array (must equal D3 cardinality in 2D) 00182 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals6, quadNodes, OPERATOR_D3), throwCounter, nException ); 00183 #endif 00184 00185 } 00186 catch (std::logic_error err) { 00187 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00188 *outStream << err.what() << '\n'; 00189 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00190 errorFlag = -1000; 00191 }; 00192 00193 // Check if number of thrown exceptions matches the one we expect 00194 if (throwCounter != nException) { 00195 errorFlag++; 00196 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00197 } 00198 00199 *outStream \ 00200 << "\n" 00201 << "===============================================================================\n"\ 00202 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\ 00203 << "===============================================================================\n"; 00204 00205 try{ 00206 std::vector<std::vector<int> > allTags = quadBasis.getAllDofTags(); 00207 00208 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again 00209 for (unsigned i = 0; i < allTags.size(); i++) { 00210 int bfOrd = quadBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]); 00211 00212 std::vector<int> myTag = quadBasis.getDofTag(bfOrd); 00213 if( !( (myTag[0] == allTags[i][0]) && 00214 (myTag[1] == allTags[i][1]) && 00215 (myTag[2] == allTags[i][2]) && 00216 (myTag[3] == allTags[i][3]) ) ) { 00217 errorFlag++; 00218 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00219 *outStream << " getDofOrdinal( {" 00220 << allTags[i][0] << ", " 00221 << allTags[i][1] << ", " 00222 << allTags[i][2] << ", " 00223 << allTags[i][3] << "}) = " << bfOrd <<" but \n"; 00224 *outStream << " getDofTag(" << bfOrd << ") = { " 00225 << myTag[0] << ", " 00226 << myTag[1] << ", " 00227 << myTag[2] << ", " 00228 << myTag[3] << "}\n"; 00229 } 00230 } 00231 00232 // Now do the same but loop over basis functions 00233 for( int bfOrd = 0; bfOrd < quadBasis.getCardinality(); bfOrd++) { 00234 std::vector<int> myTag = quadBasis.getDofTag(bfOrd); 00235 int myBfOrd = quadBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]); 00236 if( bfOrd != myBfOrd) { 00237 errorFlag++; 00238 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00239 *outStream << " getDofTag(" << bfOrd << ") = { " 00240 << myTag[0] << ", " 00241 << myTag[1] << ", " 00242 << myTag[2] << ", " 00243 << myTag[3] << "} but getDofOrdinal({" 00244 << myTag[0] << ", " 00245 << myTag[1] << ", " 00246 << myTag[2] << ", " 00247 << myTag[3] << "} ) = " << myBfOrd << "\n"; 00248 } 00249 } 00250 } 00251 catch (std::logic_error err){ 00252 *outStream << err.what() << "\n\n"; 00253 errorFlag = -1000; 00254 }; 00255 00256 *outStream \ 00257 << "\n" 00258 << "===============================================================================\n"\ 00259 << "| TEST 3: correctness of basis function values |\n"\ 00260 << "===============================================================================\n"; 00261 00262 outStream -> precision(20); 00263 00264 // VALUE: Correct basis values in (F,P) format: 00265 double basisValues[] = { 00266 1.000000000000000, 0, 0, 0, 0, 0, 0, 0, 0, -0.05333333333333334, 0, \ 00267 1.000000000000000, 0, 0, 0, 0, 0, 0, 0, 0.1066666666666667, 0, 0, \ 00268 1.000000000000000, 0, 0, 0, 0, 0, 0, -0.02666666666666666, 0, 0, 0, \ 00269 1.000000000000000, 0, 0, 0, 0, 0, 0.01333333333333333, 0, 0, 0, 0, \ 00270 1.000000000000000, 0, 0, 0, 0, 0.4266666666666667, 0, 0, 0, 0, 0, \ 00271 1.000000000000000, 0, 0, 0, 0.1422222222222222, 0, 0, 0, 0, 0, 0, \ 00272 1.000000000000000, 0, 0, -0.1066666666666667, 0, 0, 0, 0, 0, 0, 0, \ 00273 1.000000000000000, 0, -0.07111111111111112, 0, 0, 0, 0, 0, 0, 0, 0, \ 00274 1.000000000000000, 0.5688888888888890}; 00275 00276 // GRAD and D1: Correct gradients and D1 in (F,P,D) format: each row contains 6x2 values of 00277 double basisGrads[] = { 00278 -1.500000000000000, -1.500000000000000, 0.5000000000000000, 0, 0, 0, \ 00279 0, 0.5000000000000000, -0.5000000000000000, 0, 0, 0, 0, 0, 0, \ 00280 -0.5000000000000000, 0, 0, -0.08000000000000002, 0.1222222222222222, \ 00281 -0.5000000000000000, 0, 1.500000000000000, -1.500000000000000, 0, \ 00282 0.5000000000000000, 0, 0, 0.5000000000000000, 0, 0, \ 00283 -0.5000000000000000, 0, 0, 0, 0, 0, 0, 0.3999999999999999, \ 00284 -0.2444444444444444, 0, 0, 0, -0.5000000000000000, 1.500000000000000, \ 00285 1.500000000000000, -0.5000000000000000, 0, 0, 0, 0, \ 00286 0.5000000000000000, 0.5000000000000000, 0, 0, 0, 0, 0, \ 00287 -0.09999999999999998, -0.02222222222222221, 0, -0.5000000000000000, \ 00288 0, 0, 0.5000000000000000, 0, -1.500000000000000, 1.500000000000000, \ 00289 0, 0, 0, 0, -0.5000000000000000, 0, 0, 0.5000000000000000, 0, 0, \ 00290 0.02000000000000000, 0.01111111111111112, 2.000000000000000, 0, \ 00291 -2.000000000000000, 0, 0, 0, 0, 0, 0, -1.500000000000000, 0, 0, 0, \ 00292 0.5000000000000000, 0, 0, 0, -0.5000000000000000, \ 00293 -0.3199999999999999, -0.9777777777777779, 0, 0, 0, 2.000000000000000, \ 00294 0, -2.000000000000000, 0, 0, 0, 0, 1.500000000000000, 0, 0, 0, \ 00295 -0.5000000000000000, 0, 0.5000000000000000, 0, 0.5333333333333334, \ 00296 0.2666666666666666, 0, 0, 0, 0, -2.000000000000000, 0, \ 00297 2.000000000000000, 0, 0, -0.5000000000000000, 0, 0, 0, \ 00298 1.500000000000000, 0, 0, 0, 0.5000000000000000, 0.07999999999999997, \ 00299 -0.08888888888888888, 0, 2.000000000000000, 0, 0, 0, 0, 0, \ 00300 -2.000000000000000, 0, 0, 0.5000000000000000, 0, 0, 0, \ 00301 -1.500000000000000, 0, -0.5000000000000000, 0, -0.1066666666666667, \ 00302 -0.1333333333333333, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2.000000000000000, \ 00303 -2.000000000000000, 0, 0, -2.000000000000000, 2.000000000000000, 0, \ 00304 0, 0, -0.4266666666666667, 1.066666666666667}; 00305 00306 // D2: Correct multiset of second order partials in (F,P,Dk) format. D2 cardinality = 3 for 2D 00307 double basisD2[] = { 00308 1.000000000000000, 2.250000000000000, 1.000000000000000, \ 00309 1.000000000000000, -0.7500000000000000, 0, 0, 0.2500000000000000, 0, \ 00310 0, -0.7500000000000000, 1.000000000000000, 1.000000000000000, \ 00311 0.7500000000000000, 0, 0, -0.2500000000000000, 0, 0, \ 00312 -0.2500000000000000, 0, 0, 0.7500000000000000, 1.000000000000000, 0, \ 00313 0.2500000000000000, 0, 0.4800000000000000, 0.1833333333333334, \ 00314 -0.1111111111111111, 1.000000000000000, 0.7500000000000000, 0, \ 00315 1.000000000000000, -2.250000000000000, 1.000000000000000, 0, \ 00316 0.7500000000000000, 1.000000000000000, 0, -0.2500000000000000, 0, \ 00317 1.000000000000000, -0.7500000000000000, 0, 0, -0.7500000000000000, \ 00318 1.000000000000000, 0, 0.2500000000000000, 0, 0, 0.2500000000000000, \ 00319 0, 0, -0.2500000000000000, 0, 0.4800000000000000, \ 00320 -0.9166666666666666, 0.2222222222222222, 0, 0.2500000000000000, 0, 0, \ 00321 -0.7500000000000000, 1.000000000000000, 1.000000000000000, \ 00322 2.250000000000000, 1.000000000000000, 1.000000000000000, \ 00323 -0.7500000000000000, 0, 0, -0.2500000000000000, 0, 0, \ 00324 0.7500000000000000, 1.000000000000000, 1.000000000000000, \ 00325 0.7500000000000000, 0, 0, -0.2500000000000000, 0, 0, \ 00326 0.2500000000000000, 0, -0.1200000000000000, -0.08333333333333331, \ 00327 0.2222222222222222, 0, 0.7500000000000000, 1.000000000000000, 0, \ 00328 -0.2500000000000000, 0, 1.000000000000000, 0.7500000000000000, 0, \ 00329 1.000000000000000, -2.250000000000000, 1.000000000000000, 0, \ 00330 0.2500000000000000, 0, 0, 0.2500000000000000, 0, 1.000000000000000, \ 00331 -0.7500000000000000, 0, 0, -0.7500000000000000, 1.000000000000000, 0, \ 00332 -0.2500000000000000, 0, -0.1200000000000000, 0.01666666666666666, \ 00333 -0.1111111111111111, -2.000000000000000, -3.000000000000000, 0, \ 00334 -2.000000000000000, 3.000000000000000, 0, 0, -1.000000000000000, 0, \ 00335 0, 1.000000000000000, 0, -2.000000000000000, 0, 1.000000000000000, 0, \ 00336 1.000000000000000, 0, 0, 0, 1.000000000000000, 0, -1.000000000000000, \ 00337 0, 0, 0, 1.000000000000000, -0.9600000000000000, 0.7333333333333332, \ 00338 0.8888888888888890, 0, -1.000000000000000, 0, 0, 3.000000000000000, \ 00339 -2.000000000000000, 0, -3.000000000000000, -2.000000000000000, 0, \ 00340 1.000000000000000, 0, 0, 1.000000000000000, 0, 1.000000000000000, 0, \ 00341 -2.000000000000000, 0, -1.000000000000000, 0, 1.000000000000000, 0, \ 00342 0, 1.000000000000000, 0, 0, 0.6400000000000001, 1.000000000000000, \ 00343 -0.4444444444444444, 0, -1.000000000000000, 0, 0, 1.000000000000000, \ 00344 0, -2.000000000000000, -3.000000000000000, 0, -2.000000000000000, \ 00345 3.000000000000000, 0, 0, 0, 1.000000000000000, 0, -1.000000000000000, \ 00346 0, -2.000000000000000, 0, 1.000000000000000, 0, 1.000000000000000, 0, \ 00347 0, 0, 1.000000000000000, 0.2400000000000000, 0.06666666666666665, \ 00348 0.8888888888888890, 0, -3.000000000000000, -2.000000000000000, 0, \ 00349 1.000000000000000, 0, 0, -1.000000000000000, 0, 0, 3.000000000000000, \ 00350 -2.000000000000000, 0, -1.000000000000000, 0, 1.000000000000000, 0, \ 00351 0, 0, 1.000000000000000, 0, 1.000000000000000, 0, -2.000000000000000, \ 00352 1.000000000000000, 0, 0, 0.6400000000000001, -0.2000000000000001, \ 00353 0.2222222222222222, 0, 4.000000000000000, 0, 0, -4.000000000000000, \ 00354 0, 0, 4.000000000000000, 0, 0, -4.000000000000000, 0, 0, 0, \ 00355 -2.000000000000000, -2.000000000000000, 0, 0, 0, 0, \ 00356 -2.000000000000000, -2.000000000000000, 0, 0, -2.000000000000000, 0, \ 00357 -2.000000000000000, -1.280000000000000, -0.7999999999999998, \ 00358 -1.777777777777778}; 00359 00360 //D3: Correct multiset of second order partials in (F,P,Dk) format. D3 cardinality = 4 for 2D 00361 double basisD3[] = { 00362 0, -1.500000000000000, -1.500000000000000, 0, 0, -1.500000000000000, \ 00363 0.5000000000000000, 0, 0, 0.5000000000000000, 0.5000000000000000, 0, \ 00364 0, 0.5000000000000000, -1.500000000000000, 0, 0, -1.500000000000000, \ 00365 -0.5000000000000000, 0, 0, -0.5000000000000000, 0.5000000000000000, \ 00366 0, 0, 0.5000000000000000, -0.5000000000000000, 0, 0, \ 00367 -0.5000000000000000, -1.500000000000000, 0, 0, -0.5000000000000000, \ 00368 -0.5000000000000000, 0, 0, -1.100000000000000, -0.1666666666666667, \ 00369 0, 0, -1.500000000000000, -0.5000000000000000, 0, 0, \ 00370 -1.500000000000000, 1.500000000000000, 0, 0, 0.5000000000000000, \ 00371 1.500000000000000, 0, 0, 0.5000000000000000, -0.5000000000000000, 0, \ 00372 0, -1.500000000000000, 0.5000000000000000, 0, 0, -0.5000000000000000, \ 00373 1.500000000000000, 0, 0, 0.5000000000000000, 0.5000000000000000, 0, \ 00374 0, -0.5000000000000000, -0.5000000000000000, 0, 0, \ 00375 -0.5000000000000000, 0.5000000000000000, 0, 0, -1.100000000000000, \ 00376 0.8333333333333333, 0, 0, -0.5000000000000000, -0.5000000000000000, \ 00377 0, 0, -0.5000000000000000, 1.500000000000000, 0, 0, \ 00378 1.500000000000000, 1.500000000000000, 0, 0, 1.500000000000000, \ 00379 -0.5000000000000000, 0, 0, -0.5000000000000000, 0.5000000000000000, \ 00380 0, 0, 0.5000000000000000, 1.500000000000000, 0, 0, 1.500000000000000, \ 00381 0.5000000000000000, 0, 0, 0.5000000000000000, -0.5000000000000000, 0, \ 00382 0, 0.5000000000000000, 0.5000000000000000, 0, 0, \ 00383 -0.09999999999999998, 0.8333333333333333, 0, 0, -0.5000000000000000, \ 00384 -1.500000000000000, 0, 0, -0.5000000000000000, 0.5000000000000000, 0, \ 00385 0, 1.500000000000000, 0.5000000000000000, 0, 0, 1.500000000000000, \ 00386 -1.500000000000000, 0, 0, -0.5000000000000000, -0.5000000000000000, \ 00387 0, 0, 0.5000000000000000, 0.5000000000000000, 0, 0, \ 00388 1.500000000000000, -0.5000000000000000, 0, 0, 0.5000000000000000, \ 00389 -1.500000000000000, 0, 0, 0.5000000000000000, -0.5000000000000000, 0, \ 00390 0, -0.09999999999999998, -0.1666666666666667, 0, 0, \ 00391 3.000000000000000, 2.000000000000000, 0, 0, 3.000000000000000, \ 00392 -2.000000000000000, 0, 0, -1.000000000000000, -2.000000000000000, 0, \ 00393 0, -1.000000000000000, 2.000000000000000, 0, 0, 3.000000000000000, 0, \ 00394 0, 0, 1.000000000000000, -2.000000000000000, 0, 0, \ 00395 -1.000000000000000, 0, 0, 0, 1.000000000000000, 2.000000000000000, 0, \ 00396 0, 1.000000000000000, 0, 0, 0, 2.200000000000000, \ 00397 -0.6666666666666665, 0, 0, 2.000000000000000, 1.000000000000000, 0, \ 00398 0, 2.000000000000000, -3.000000000000000, 0, 0, -2.000000000000000, \ 00399 -3.000000000000000, 0, 0, -2.000000000000000, 1.000000000000000, 0, \ 00400 0, 2.000000000000000, -1.000000000000000, 0, 0, 0, \ 00401 -3.000000000000000, 0, 0, -2.000000000000000, -1.000000000000000, 0, \ 00402 0, 0, 1.000000000000000, 0, 0, 0, -1.000000000000000, 0, 0, \ 00403 1.200000000000000, -1.666666666666667, 0, 0, 1.000000000000000, \ 00404 2.000000000000000, 0, 0, 1.000000000000000, -2.000000000000000, 0, 0, \ 00405 -3.000000000000000, -2.000000000000000, 0, 0, -3.000000000000000, \ 00406 2.000000000000000, 0, 0, 1.000000000000000, 0, 0, 0, \ 00407 -1.000000000000000, -2.000000000000000, 0, 0, -3.000000000000000, 0, \ 00408 0, 0, -1.000000000000000, 2.000000000000000, 0, 0, \ 00409 -1.000000000000000, 0, 0, 0, 0.2000000000000000, -0.6666666666666665, \ 00410 0, 0, 2.000000000000000, 3.000000000000000, 0, 0, 2.000000000000000, \ 00411 -1.000000000000000, 0, 0, -2.000000000000000, -1.000000000000000, 0, \ 00412 0, -2.000000000000000, 3.000000000000000, 0, 0, 2.000000000000000, \ 00413 1.000000000000000, 0, 0, 0, -1.000000000000000, 0, 0, \ 00414 -2.000000000000000, 1.000000000000000, 0, 0, 0, 3.000000000000000, 0, \ 00415 0, 0, 1.000000000000000, 0, 0, 1.200000000000000, 0.3333333333333334, \ 00416 0, 0, -4.000000000000000, -4.000000000000000, 0, 0, \ 00417 -4.000000000000000, 4.000000000000000, 0, 0, 4.000000000000000, \ 00418 4.000000000000000, 0, 0, 4.000000000000000, -4.000000000000000, 0, 0, \ 00419 -4.000000000000000, 0, 0, 0, 0, 4.000000000000000, 0, 0, \ 00420 4.000000000000000, 0, 0, 0, 0, -4.000000000000000, 0, 0, 0, 0, 0, 0, \ 00421 -2.400000000000000, 1.333333333333333, 0}; 00422 //D4: Correct multiset of second order partials in (F,P,Dk) format. D4 cardinality = 5 for 2D 00423 double basisD4[] = { 00424 0, 0, 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00425 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00426 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00427 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00428 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00429 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00430 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00431 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00432 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00433 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00434 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00435 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00436 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00437 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00438 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00439 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00440 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00441 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00442 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00443 1.000000000000000, 0, 0, 0, 0, 1.000000000000000, 0, 0, 0, 0, \ 00444 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00445 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00446 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00447 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00448 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00449 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00450 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00451 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00452 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00453 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00454 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00455 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00456 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00457 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00458 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00459 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00460 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00461 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00462 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00463 -2.000000000000000, 0, 0, 0, 0, -2.000000000000000, 0, 0, 0, 0, \ 00464 4.000000000000000, 0, 0, 0, 0, 4.000000000000000, 0, 0, 0, 0, \ 00465 4.000000000000000, 0, 0, 0, 0, 4.000000000000000, 0, 0, 0, 0, \ 00466 4.000000000000000, 0, 0, 0, 0, 4.000000000000000, 0, 0, 0, 0, \ 00467 4.000000000000000, 0, 0, 0, 0, 4.000000000000000, 0, 0, 0, 0, \ 00468 4.000000000000000, 0, 0, 0, 0, 4.000000000000000, 0, 0}; 00469 00470 try{ 00471 00472 // Dimensions for the output arrays: 00473 int numFields = quadBasis.getCardinality(); 00474 int numPoints = quadNodes.dimension(0); 00475 int spaceDim = quadBasis.getBaseCellTopology().getDimension(); 00476 00477 // Generic array for values, grads, curls, etc. that will be properly sized before each call 00478 FieldContainer<double> vals; 00479 00480 // Check VALUE of basis functions: resize vals to rank-2 container: 00481 vals.resize(numFields, numPoints); 00482 quadBasis.getValues(vals, quadNodes, OPERATOR_VALUE); 00483 for (int i = 0; i < numFields; i++) { 00484 for (int j = 0; j < numPoints; j++) { 00485 00486 // Compute offset for (F,P) container 00487 int l = j + i * numPoints; 00488 if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) { 00489 errorFlag++; 00490 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00491 00492 // Output the multi-index of the value where the error is: 00493 *outStream << " At multi-index { "; 00494 *outStream << i << " ";*outStream << j << " "; 00495 *outStream << "} computed value: " << vals(i,j) 00496 << " but reference value: " << basisValues[l] << "\n"; 00497 } 00498 } 00499 } 00500 00501 // Check GRAD of basis function: resize vals to rank-3 container 00502 vals.resize(numFields, numPoints, spaceDim); 00503 quadBasis.getValues(vals, quadNodes, OPERATOR_GRAD); 00504 for (int i = 0; i < numFields; i++) { 00505 for (int j = 0; j < numPoints; j++) { 00506 for (int k = 0; k < spaceDim; k++) { 00507 00508 // basisGrads is (F,P,D), compute offset: 00509 int l = k + j * spaceDim + i * spaceDim * numPoints; 00510 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00511 errorFlag++; 00512 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00513 00514 // Output the multi-index of the value where the error is: 00515 *outStream << " At multi-index { "; 00516 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00517 *outStream << "} computed grad component: " << vals(i,j,k) 00518 << " but reference grad component: " << basisGrads[l] << "\n"; 00519 } 00520 } 00521 } 00522 } 00523 00524 00525 // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD) 00526 quadBasis.getValues(vals, quadNodes, OPERATOR_D1); 00527 for (int i = 0; i < numFields; i++) { 00528 for (int j = 0; j < numPoints; j++) { 00529 for (int k = 0; k < spaceDim; k++) { 00530 00531 // basisGrads is (F,P,D), compute offset: 00532 int l = k + j * spaceDim + i * spaceDim * numPoints; 00533 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) { 00534 errorFlag++; 00535 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00536 00537 // Output the multi-index of the value where the error is: 00538 *outStream << " At multi-index { "; 00539 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00540 *outStream << "} computed D1 component: " << vals(i,j,k) 00541 << " but reference D1 component: " << basisGrads[l] << "\n"; 00542 } 00543 } 00544 } 00545 } 00546 00547 00548 // Check CURL of basis function: resize vals just for illustration! 00549 vals.resize(numFields, numPoints, spaceDim); 00550 quadBasis.getValues(vals, quadNodes, OPERATOR_CURL); 00551 for (int i = 0; i < numFields; i++) { 00552 for (int j = 0; j < numPoints; j++) { 00553 // We will use "rotated" basisGrads to check CURL: get offsets to extract (u_y, -u_x) 00554 int curl_0 = 1 + j * spaceDim + i * spaceDim * numPoints; // position of y-derivative 00555 int curl_1 = 0 + j * spaceDim + i * spaceDim * numPoints; // position of x-derivative 00556 00557 double curl_value_0 = basisGrads[curl_0]; 00558 double curl_value_1 =-basisGrads[curl_1]; 00559 if (std::abs(vals(i,j,0) - curl_value_0) > INTREPID_TOL) { 00560 errorFlag++; 00561 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00562 // Output the multi-index of the value where the error is: 00563 *outStream << " At multi-index { "; 00564 *outStream << i << " ";*outStream << j << " ";*outStream << 0 << " "; 00565 *outStream << "} computed curl component: " << vals(i,j,0) 00566 << " but reference curl component: " << curl_value_0 << "\n"; 00567 } 00568 if (std::abs(vals(i,j,1) - curl_value_1) > INTREPID_TOL) { 00569 errorFlag++; 00570 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00571 // Output the multi-index of the value where the error is: 00572 *outStream << " At multi-index { "; 00573 *outStream << i << " ";*outStream << j << " ";*outStream << 1 << " "; 00574 *outStream << "} computed curl component: " << vals(i,j,1) 00575 << " but reference curl component: " << curl_value_1 << "\n"; 00576 } 00577 } 00578 } 00579 00580 // Check D2 of basis function 00581 int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim); 00582 vals.resize(numFields, numPoints, D2cardinality); 00583 quadBasis.getValues(vals, quadNodes, OPERATOR_D2); 00584 for (int i = 0; i < numFields; i++) { 00585 for (int j = 0; j < numPoints; j++) { 00586 for (int k = 0; k < D2cardinality; k++) { 00587 00588 // basisD2 is (F,P,Dk), compute offset: 00589 int l = k + j * D2cardinality + i * D2cardinality * numPoints; 00590 if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) { 00591 errorFlag++; 00592 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00593 00594 // Output the multi-index of the value where the error is: 00595 *outStream << " At multi-index { "; 00596 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00597 *outStream << "} computed D2 component: " << vals(i,j,k) 00598 << " but reference D2 component: " << basisD2[l] << "\n"; 00599 } 00600 } 00601 } 00602 } 00603 00604 00605 // Check D3 of basis function 00606 int D3cardinality = Intrepid::getDkCardinality(OPERATOR_D3, spaceDim); 00607 vals.resize(numFields, numPoints, D3cardinality); 00608 quadBasis.getValues(vals, quadNodes, OPERATOR_D3); 00609 for (int i = 0; i < numFields; i++) { 00610 for (int j = 0; j < numPoints; j++) { 00611 for (int k = 0; k < D3cardinality; k++) { 00612 00613 // basisD3 is (F,P,Dk), compute offset: 00614 int l = k + j * D3cardinality + i * D3cardinality * numPoints; 00615 if (std::abs(vals(i,j,k) - basisD3[l]) > INTREPID_TOL) { 00616 errorFlag++; 00617 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00618 00619 // Output the multi-index of the value where the error is: 00620 *outStream << " At multi-index { "; 00621 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00622 *outStream << "} computed D3 component: " << vals(i,j,k) 00623 << " but reference D3 component: " << basisD2[l] << "\n"; 00624 } 00625 } 00626 } 00627 } 00628 00629 00630 // Check D4 of basis function 00631 int D4cardinality = Intrepid::getDkCardinality(OPERATOR_D4, spaceDim); 00632 vals.resize(numFields, numPoints, D4cardinality); 00633 quadBasis.getValues(vals, quadNodes, OPERATOR_D4); 00634 for (int i = 0; i < numFields; i++) { 00635 for (int j = 0; j < numPoints; j++) { 00636 for (int k = 0; k < D4cardinality; k++) { 00637 00638 // basisD4 is (F,P,Dk), compute offset: 00639 int l = k + j * D4cardinality + i * D4cardinality * numPoints; 00640 if (std::abs(vals(i,j,k) - basisD4[l]) > INTREPID_TOL) { 00641 errorFlag++; 00642 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00643 00644 // Output the multi-index of the value where the error is: 00645 *outStream << " At multi-index { "; 00646 *outStream << i << " ";*outStream << j << " ";*outStream << k << " "; 00647 *outStream << "} computed D4 component: " << vals(i,j,k) 00648 << " but reference D4 component: " << basisD2[l] << "\n"; 00649 } 00650 } 00651 } 00652 } 00653 00654 00655 // Check all higher derivatives - must be zero. 00656 for(EOperator op = OPERATOR_D5; op < OPERATOR_MAX; op++) { 00657 00658 // The last dimension is the number of kth derivatives and needs to be resized for every Dk 00659 int DkCardin = Intrepid::getDkCardinality(op, spaceDim); 00660 vals.resize(numFields, numPoints, DkCardin); 00661 00662 quadBasis.getValues(vals, quadNodes, op); 00663 for (int i = 0; i < vals.size(); i++) { 00664 if (std::abs(vals[i]) > INTREPID_TOL) { 00665 errorFlag++; 00666 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00667 00668 // Get the multi-index of the value where the error is and the operator order 00669 std::vector<int> myIndex; 00670 vals.getMultiIndex(myIndex,i); 00671 int ord = Intrepid::getOperatorOrder(op); 00672 *outStream << " At multi-index { "; 00673 for(int j = 0; j < vals.rank(); j++) { 00674 *outStream << myIndex[j] << " "; 00675 } 00676 *outStream << "} computed D"<< ord <<" component: " << vals[i] 00677 << " but reference D" << ord << " component: 0 \n"; 00678 } 00679 } 00680 } 00681 } 00682 // Catch unexpected errors 00683 catch (std::logic_error err) { 00684 *outStream << err.what() << "\n\n"; 00685 errorFlag = -1000; 00686 }; 00687 00688 00689 *outStream \ 00690 << "\n" 00691 << "===============================================================================\n"\ 00692 << "| TEST 4: correctness of DoF locations |\n"\ 00693 << "===============================================================================\n"; 00694 00695 try{ 00696 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis = 00697 Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM<double, FieldContainer<double> >); 00698 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface = 00699 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis); 00700 00701 FieldContainer<double> cvals; 00702 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality()); 00703 00704 // Check exceptions. 00705 #ifdef HAVE_INTREPID_DEBUG 00706 cvals.resize(1,2,3); 00707 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00708 cvals.resize(8,2); 00709 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00710 cvals.resize(9,1); 00711 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); 00712 #endif 00713 cvals.resize(9,2); 00714 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--; 00715 // Check if number of thrown exceptions matches the one we expect 00716 if (throwCounter != nException) { 00717 errorFlag++; 00718 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00719 } 00720 00721 // Check mathematical correctness. 00722 basis->getValues(bvals, cvals, OPERATOR_VALUE); 00723 char buffer[120]; 00724 for (int i=0; i<bvals.dimension(0); i++) { 00725 for (int j=0; j<bvals.dimension(1); j++) { 00726 if ((i != j) && (std::abs(bvals(i,j) - 0.0) > INTREPID_TOL)) { 00727 errorFlag++; 00728 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), bvals(i,j), 0.0); 00729 *outStream << buffer; 00730 } 00731 else if ((i == j) && (std::abs(bvals(i,j) - 1.0) > INTREPID_TOL)) { 00732 errorFlag++; 00733 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), bvals(i,j), 1.0); 00734 *outStream << buffer; 00735 } 00736 } 00737 } 00738 00739 } 00740 catch (std::logic_error err){ 00741 *outStream << err.what() << "\n\n"; 00742 errorFlag = -1000; 00743 }; 00744 00745 if (errorFlag != 0) 00746 std::cout << "End Result: TEST FAILED\n"; 00747 else 00748 std::cout << "End Result: TEST PASSED\n"; 00749 00750 // reset format state of std::cout 00751 std::cout.copyfmt(oldFormatState); 00752 00753 return errorFlag; 00754 }
1.7.4