|
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) or 00025 // Denis Ridzal (dridzal@sandia.gov). 00026 // 00027 // ************************************************************************ 00028 // @HEADER 00029 00035 #include "Intrepid_FunctionSpaceTools.hpp" 00036 #include "Intrepid_FieldContainer.hpp" 00037 #include "Intrepid_CellTools.hpp" 00038 #include "Intrepid_HGRAD_TET_C1_FEM.hpp" 00039 #include "Intrepid_DefaultCubatureFactory.hpp" 00040 #include "Intrepid_Utils.hpp" 00041 #include "Teuchos_oblackholestream.hpp" 00042 #include "Teuchos_RCP.hpp" 00043 #include "Teuchos_GlobalMPISession.hpp" 00044 00045 using namespace std; 00046 using namespace Intrepid; 00047 00048 #define INTREPID_TEST_COMMAND( S ) \ 00049 { \ 00050 try { \ 00051 S ; \ 00052 } \ 00053 catch (std::logic_error err) { \ 00054 *outStream << "Expected Error ----------------------------------------------------------------\n"; \ 00055 *outStream << err.what() << '\n'; \ 00056 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \ 00057 }; \ 00058 } 00059 00060 00061 int main(int argc, char *argv[]) { 00062 00063 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00064 00065 // This little trick lets us print to std::cout only if 00066 // a (dummy) command-line argument is provided. 00067 int iprint = argc - 1; 00068 Teuchos::RCP<std::ostream> outStream; 00069 Teuchos::oblackholestream bhs; // outputs nothing 00070 if (iprint > 0) 00071 outStream = Teuchos::rcp(&std::cout, false); 00072 else 00073 outStream = Teuchos::rcp(&bhs, false); 00074 00075 // Save the format state of the original std::cout. 00076 Teuchos::oblackholestream oldFormatState; 00077 oldFormatState.copyfmt(std::cout); 00078 00079 *outStream \ 00080 << "===============================================================================\n" \ 00081 << "| |\n" \ 00082 << "| Unit Test (FunctionSpaceTools) |\n" \ 00083 << "| |\n" \ 00084 << "| 1) basic operator transformations and integration in HGRAD |\n" \ 00085 << "| |\n" \ 00086 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \ 00087 << "| Denis Ridzal (dridzal@sandia.gov). |\n" \ 00088 << "| |\n" \ 00089 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00090 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00091 << "| |\n" \ 00092 << "===============================================================================\n"; 00093 00094 00095 int errorFlag = 0; 00096 #ifdef HAVE_INTREPID_DEBUG 00097 int beginThrowNumber = TestForException_getThrowNumber(); 00098 int endThrowNumber = beginThrowNumber + 28; 00099 #endif 00100 00101 typedef FunctionSpaceTools fst; 00102 00103 *outStream \ 00104 << "\n" 00105 << "===============================================================================\n"\ 00106 << "| TEST 1: exceptions |\n"\ 00107 << "===============================================================================\n"; 00108 00109 try{ 00110 #ifdef HAVE_INTREPID_DEBUG 00111 FieldContainer<double> a_2(2); 00112 FieldContainer<double> a_2_2(2, 2); 00113 FieldContainer<double> a_2_3(2, 3); 00114 FieldContainer<double> a_3_2(3, 2); 00115 FieldContainer<double> a_2_2_3(2, 2, 3); 00116 FieldContainer<double> a_2_2_3_3(2, 2, 3, 3); 00117 FieldContainer<double> a_2_2_2(2, 2, 2); 00118 FieldContainer<double> a_2_2_2_3_3(2, 2, 2, 3, 3); 00119 FieldContainer<double> a_2_2_2_2_2(2, 2, 2, 2, 2); 00120 FieldContainer<double> a_2_2_2_2(2, 2, 2, 2); 00121 FieldContainer<double> a_3_2_2_2(3, 2, 2, 2); 00122 FieldContainer<double> a_2_3_2_2(2, 3, 2, 2); 00123 FieldContainer<double> a_2_2_3_2(2, 2, 3, 2); 00124 FieldContainer<double> a_2_2_2_3(2, 2, 2, 3); 00125 00126 *outStream << "\n >>>>> TESTING computeCellMeasure:\n"; 00127 INTREPID_TEST_COMMAND( fst::computeCellMeasure<double>(a_2_2, a_2, a_2) ); 00128 INTREPID_TEST_COMMAND( fst::computeCellMeasure<double>(a_2_2, a_2_2, a_2) ); 00129 00130 *outStream << "\n >>>>> TESTING computeFaceMeasure:\n"; 00131 INTREPID_TEST_COMMAND( fst::computeFaceMeasure<double>(a_2_2, a_2, a_2, 0, shards::getCellTopologyData< shards::Tetrahedron<> >()) ); 00132 INTREPID_TEST_COMMAND( fst::computeFaceMeasure<double>(a_2_2, a_2_2_3_3, a_2, 0, shards::getCellTopologyData< shards::Tetrahedron<> >()) ); 00133 00134 *outStream << "\n >>>>> TESTING computeEdgeMeasure:\n"; 00135 INTREPID_TEST_COMMAND( fst::computeEdgeMeasure<double>(a_2_2, a_2, a_2, 0, shards::getCellTopologyData< shards::Triangle<> >()) ); 00136 INTREPID_TEST_COMMAND( fst::computeEdgeMeasure<double>(a_2_2, a_2_2_2_2, a_2, 0, shards::getCellTopologyData< shards::Triangle<> >()) ); 00137 00138 *outStream << "\n >>>>> TESTING integrate:\n"; 00139 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2_2_2_2, a_2_2_2, a_2_2_2, COMP_CPP) ); 00140 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2, a_2_2, a_2_2, COMP_CPP) ); 00141 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2, a_2_2_3, a_2_2_3, COMP_CPP) ); 00142 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2, a_2_2_3_3, a_2_2_3_3, COMP_CPP) ); 00143 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2_2, a_2_2, a_2_2_2, COMP_CPP) ); 00144 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2_2, a_2_2_3, a_2_2_2_3, COMP_CPP) ); 00145 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2_2, a_2_2_3_3, a_2_2_2_3_3, COMP_CPP) ); 00146 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2_2_2, a_2_2_2, a_2_2_2, COMP_CPP) ); 00147 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2_2_2, a_2_2_2_3, a_2_2_2_3, COMP_CPP) ); 00148 INTREPID_TEST_COMMAND( fst::integrate<double>(a_2_2_2, a_2_2_2_3_3, a_2_2_2_3_3, COMP_CPP) ); 00149 00150 *outStream << "\n >>>>> TESTING operatorIntegral:\n"; 00151 INTREPID_TEST_COMMAND( fst::operatorIntegral<double>(a_2_2_2, a_2_2, a_2_2_2, COMP_CPP) ); 00152 INTREPID_TEST_COMMAND( fst::operatorIntegral<double>(a_2_2_2, a_2_2_2, a_2_2_2, COMP_CPP) ); 00153 INTREPID_TEST_COMMAND( fst::operatorIntegral<double>(a_2_2_2, a_2_2_2_3, a_2_2_2_3, COMP_CPP) ); 00154 INTREPID_TEST_COMMAND( fst::operatorIntegral<double>(a_2_2_2, a_2_2_2_3_3, a_2_2_2_3_3, COMP_CPP) ); 00155 00156 *outStream << "\n >>>>> TESTING functionalIntegral:\n"; 00157 INTREPID_TEST_COMMAND( fst::functionalIntegral<double>(a_2_2, a_2_2_2_3_3, a_2_2_2, COMP_CPP) ); 00158 INTREPID_TEST_COMMAND( fst::functionalIntegral<double>(a_2_2, a_2_2, a_2_2_2, COMP_CPP) ); 00159 INTREPID_TEST_COMMAND( fst::functionalIntegral<double>(a_2_2, a_2_2_3, a_2_2_2_3, COMP_CPP) ); 00160 INTREPID_TEST_COMMAND( fst::functionalIntegral<double>(a_2_2, a_2_2_3_3, a_2_2_2_3_3, COMP_CPP) ); 00161 00162 *outStream << "\n >>>>> TESTING dataIntegral:\n"; 00163 INTREPID_TEST_COMMAND( fst::dataIntegral<double>(a_2, a_2, a_2_2_2, COMP_CPP) ); 00164 INTREPID_TEST_COMMAND( fst::dataIntegral<double>(a_2, a_2_2, a_2_2, COMP_CPP) ); 00165 INTREPID_TEST_COMMAND( fst::dataIntegral<double>(a_2, a_2_2_3, a_2_2_3, COMP_CPP) ); 00166 INTREPID_TEST_COMMAND( fst::dataIntegral<double>(a_2, a_2_2_3_3, a_2_2_3_3, COMP_CPP) ); 00167 00168 *outStream << "\n >>>>> TESTING applyLeftFieldSigns:\n"; 00169 INTREPID_TEST_COMMAND( fst::applyLeftFieldSigns<double>(a_2, a_2) ); 00170 INTREPID_TEST_COMMAND( fst::applyLeftFieldSigns<double>(a_2_2_2, a_2) ); 00171 INTREPID_TEST_COMMAND( fst::applyLeftFieldSigns<double>(a_2_2_2, a_3_2) ); 00172 INTREPID_TEST_COMMAND( fst::applyLeftFieldSigns<double>(a_2_2_2, a_2_3) ); 00173 INTREPID_TEST_COMMAND( fst::applyLeftFieldSigns<double>(a_2_2_2, a_2_2) ); 00174 00175 *outStream << "\n >>>>> TESTING applyRightFieldSigns:\n"; 00176 INTREPID_TEST_COMMAND( fst::applyRightFieldSigns<double>(a_2, a_2) ); 00177 INTREPID_TEST_COMMAND( fst::applyRightFieldSigns<double>(a_2_2_2, a_2) ); 00178 INTREPID_TEST_COMMAND( fst::applyRightFieldSigns<double>(a_2_2_2, a_3_2) ); 00179 INTREPID_TEST_COMMAND( fst::applyRightFieldSigns<double>(a_2_2_2, a_2_3) ); 00180 INTREPID_TEST_COMMAND( fst::applyRightFieldSigns<double>(a_2_2_2, a_2_2) ); 00181 00182 *outStream << "\n >>>>> TESTING applyFieldSigns:\n"; 00183 INTREPID_TEST_COMMAND( fst::applyFieldSigns<double>(a_2, a_2) ); 00184 INTREPID_TEST_COMMAND( fst::applyFieldSigns<double>(a_2_2, a_2) ); 00185 INTREPID_TEST_COMMAND( fst::applyFieldSigns<double>(a_2_2, a_3_2) ); 00186 INTREPID_TEST_COMMAND( fst::applyFieldSigns<double>(a_2_2, a_2_3) ); 00187 INTREPID_TEST_COMMAND( fst::applyFieldSigns<double>(a_2_2_2_3_3, a_2_2) ); 00188 00189 *outStream << "\n >>>>> TESTING evaluate:\n"; 00190 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2, a_2, a_2_2) ); 00191 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2, a_2, a_2_2_2_3_3) ); 00192 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2, a_2_2, a_2_2_2_3_3) ); 00193 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2_2_3_3, a_3_2, a_2_2_2_3_3) ); 00194 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2_2_3_3, a_2_3, a_2_2_2_3_3) ); 00195 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_3_2_2_2, a_2_2, a_2_2_2_2_2) ); 00196 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2_3_2_2, a_2_2, a_2_2_2_2_2) ); 00197 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2_2_3_2, a_2_2, a_2_2_2_2_2) ); 00198 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2_2_2_3, a_2_2, a_2_2_2_2_2) ); 00199 INTREPID_TEST_COMMAND( fst::evaluate<double>(a_2_2_2_2, a_2_2, a_2_2_2_2_2) ); 00200 #endif 00201 } 00202 catch (std::logic_error err) { 00203 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00204 *outStream << err.what() << '\n'; 00205 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00206 errorFlag = -1000; 00207 }; 00208 00209 #ifdef HAVE_INTREPID_DEBUG 00210 if (TestForException_getThrowNumber() != endThrowNumber) 00211 errorFlag++; 00212 #endif 00213 00214 *outStream \ 00215 << "\n" 00216 << "===============================================================================\n"\ 00217 << "| TEST 2: correctness of math operations |\n"\ 00218 << "===============================================================================\n"; 00219 00220 outStream->precision(20); 00221 00222 try { 00223 shards::CellTopology cellType = shards::getCellTopologyData< shards::Tetrahedron<> >(); // cell type: tet 00224 00225 /* Related to cubature. */ 00226 DefaultCubatureFactory<double> cubFactory; // create cubature factory 00227 int cubDegree = 2; // cubature degree 00228 Teuchos::RCP<Cubature<double> > myCub = cubFactory.create(cellType, cubDegree); // create default cubature 00229 int spaceDim = myCub->getDimension(); // get spatial dimension 00230 int numCubPoints = myCub->getNumPoints(); // get number of cubature points 00231 00232 /* Related to basis. */ 00233 Basis_HGRAD_TET_C1_FEM<double, FieldContainer<double> > tetBasis; // create tet basis 00234 int numFields = tetBasis.getCardinality(); // get basis cardinality 00235 00236 /* Cell geometries. */ 00237 int numCells = 4; 00238 int numNodes = 4; 00239 int numCellData = numCells*numNodes*spaceDim; 00240 double tetnodes[] = { 00241 // tet 0 00242 0.0, 0.0, 0.0, 00243 1.0, 0.0, 0.0, 00244 0.0, 1.0, 0.0, 00245 0.0, 0.0, 1.0, 00246 // tet 1 00247 4.0, 5.0, 1.0, 00248 -6.0, 2.0, 0.0, 00249 4.0, -3.0, -1.0, 00250 0.0, 2.0, 5.0, 00251 // tet 2 00252 -6.0, -3.0, 1.0, 00253 9.0, 2.0, 1.0, 00254 8.9, 2.1, 0.9, 00255 8.9, 2.1, 1.1, 00256 // tet 3 00257 -6.0, -3.0, 1.0, 00258 12.0, 3.0, 1.0, 00259 2.9, 0.1, 0.9, 00260 2.9, 0.1, 1.1 00261 }; 00262 00263 /* Computational arrays. */ 00264 FieldContainer<double> cub_points(numCubPoints, spaceDim); 00265 FieldContainer<double> cub_weights(numCubPoints); 00266 FieldContainer<double> cell_nodes(numCells, numNodes, spaceDim); 00267 FieldContainer<double> jacobian(numCells, numCubPoints, spaceDim, spaceDim); 00268 FieldContainer<double> jacobian_inv(numCells, numCubPoints, spaceDim, spaceDim); 00269 FieldContainer<double> jacobian_det(numCells, numCubPoints); 00270 FieldContainer<double> weighted_measure(numCells, numCubPoints); 00271 00272 FieldContainer<double> grad_of_basis_at_cub_points(numFields, numCubPoints, spaceDim); 00273 FieldContainer<double> transformed_grad_of_basis_at_cub_points(numCells, numFields, numCubPoints, spaceDim); 00274 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points(numCells, numFields, numCubPoints, spaceDim); 00275 FieldContainer<double> stiffness_matrices(numCells, numFields, numFields); 00276 00277 FieldContainer<double> value_of_basis_at_cub_points(numFields, numCubPoints); 00278 FieldContainer<double> transformed_value_of_basis_at_cub_points(numCells, numFields, numCubPoints); 00279 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points(numCells, numFields, numCubPoints); 00280 FieldContainer<double> mass_matrices(numCells, numFields, numFields); 00281 00282 /******************* START COMPUTATION ***********************/ 00283 00284 // get cubature points and weights 00285 myCub->getCubature(cub_points, cub_weights); 00286 00287 // fill cell vertex array 00288 cell_nodes.setValues(tetnodes, numCellData); 00289 00290 // compute geometric cell information 00291 CellTools<double>::setJacobian(jacobian, cub_points, cell_nodes, cellType); 00292 CellTools<double>::setJacobianInv(jacobian_inv, jacobian); 00293 CellTools<double>::setJacobianDet(jacobian_det, jacobian); 00294 00295 // compute weighted measure 00296 fst::computeCellMeasure<double>(weighted_measure, jacobian_det, cub_weights); 00297 00298 // Computing stiffness matrices: 00299 // tabulate gradients of basis functions at (reference) cubature points 00300 tetBasis.getValues(grad_of_basis_at_cub_points, cub_points, OPERATOR_GRAD); 00301 00302 // transform gradients of basis functions 00303 fst::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points, 00304 jacobian_inv, 00305 grad_of_basis_at_cub_points); 00306 00307 // multiply with weighted measure 00308 fst::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points, 00309 weighted_measure, 00310 transformed_grad_of_basis_at_cub_points); 00311 00312 // compute stiffness matrices 00313 fst::integrate<double>(stiffness_matrices, 00314 transformed_grad_of_basis_at_cub_points, 00315 weighted_transformed_grad_of_basis_at_cub_points, 00316 COMP_CPP); 00317 00318 00319 // Computing mass matrices: 00320 // tabulate values of basis functions at (reference) cubature points 00321 tetBasis.getValues(value_of_basis_at_cub_points, cub_points, OPERATOR_VALUE); 00322 00323 // transform values of basis functions 00324 fst::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points, 00325 value_of_basis_at_cub_points); 00326 00327 // multiply with weighted measure 00328 fst::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points, 00329 weighted_measure, 00330 transformed_value_of_basis_at_cub_points); 00331 00332 // compute mass matrices 00333 fst::integrate<double>(mass_matrices, 00334 transformed_value_of_basis_at_cub_points, 00335 weighted_transformed_value_of_basis_at_cub_points, 00336 COMP_CPP); 00337 00338 /******************* STOP COMPUTATION ***********************/ 00339 00340 00341 /******************* START COMPARISON ***********************/ 00342 string basedir = "./testdata"; 00343 for (int cell_id = 0; cell_id < numCells; cell_id++) { 00344 00345 stringstream namestream; 00346 string filename; 00347 namestream << basedir << "/mass_TET_FEM_P1" << "_" << "0" << cell_id+1 << ".dat"; 00348 namestream >> filename; 00349 00350 ifstream massfile(&filename[0]); 00351 if (massfile.is_open()) { 00352 if (compareToAnalytic<double>(&mass_matrices(cell_id, 0, 0), massfile, 1e-10, 0) > 0) 00353 errorFlag++; 00354 massfile.close(); 00355 } 00356 else { 00357 errorFlag = -1; 00358 std::cout << "End Result: TEST FAILED\n"; 00359 return errorFlag; 00360 } 00361 00362 namestream.clear(); 00363 namestream << basedir << "/stiff_TET_FEM_P1" << "_" << "0" << cell_id+1 << ".dat"; 00364 namestream >> filename; 00365 00366 ifstream stifffile(&filename[0]); 00367 if (stifffile.is_open()) 00368 { 00369 if (compareToAnalytic<double>(&stiffness_matrices(cell_id, 0, 0), stifffile, 1e-10, 0) > 0) 00370 errorFlag++; 00371 stifffile.close(); 00372 } 00373 else { 00374 errorFlag = -1; 00375 std::cout << "End Result: TEST FAILED\n"; 00376 return errorFlag; 00377 } 00378 00379 } 00380 00381 /******************* STOP COMPARISON ***********************/ 00382 00383 *outStream << "\n"; 00384 } 00385 catch (std::logic_error err) { 00386 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n"; 00387 *outStream << err.what() << '\n'; 00388 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; 00389 errorFlag = -1000; 00390 }; 00391 00392 00393 if (errorFlag != 0) 00394 std::cout << "End Result: TEST FAILED\n"; 00395 else 00396 std::cout << "End Result: TEST PASSED\n"; 00397 00398 // reset format state of std::cout 00399 std::cout.copyfmt(oldFormatState); 00400 00401 return errorFlag; 00402 }
1.7.4