|
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_LINE_Cn_FEM.hpp" 00038 #include "Intrepid_DefaultCubatureFactory.hpp" 00039 #include "Intrepid_RealSpaceTools.hpp" 00040 #include "Intrepid_ArrayTools.hpp" 00041 #include "Intrepid_FunctionSpaceTools.hpp" 00042 #include "Intrepid_CellTools.hpp" 00043 #include "Intrepid_PointTools.hpp" 00044 #include "Teuchos_oblackholestream.hpp" 00045 #include "Teuchos_RCP.hpp" 00046 #include "Teuchos_GlobalMPISession.hpp" 00047 #include "Teuchos_SerialDenseMatrix.hpp" 00048 #include "Teuchos_SerialDenseVector.hpp" 00049 #include "Teuchos_LAPACK.hpp" 00050 00051 using namespace std; 00052 using namespace Intrepid; 00053 00054 void rhsFunc(FieldContainer<double> &, const FieldContainer<double> &, int); 00055 void neumann(FieldContainer<double> &, const FieldContainer<double> &, const FieldContainer<double> &, int); 00056 void u_exact(FieldContainer<double> &, const FieldContainer<double> &, int); 00057 00059 void rhsFunc(FieldContainer<double> & result, const FieldContainer<double> & points, int degree) { 00060 if (degree == 0) { 00061 for (int cell=0; cell<result.dimension(0); cell++) { 00062 for (int pt=0; pt<result.dimension(1); pt++) { 00063 result(cell,pt) = 1.0; 00064 00065 } 00066 } 00067 } 00068 else if (degree == 1) { 00069 for (int cell=0; cell<result.dimension(0); cell++) { 00070 for (int pt=0; pt<result.dimension(1); pt++) { 00071 result(cell,pt) = points(cell,pt,0); 00072 } 00073 } 00074 } 00075 else { 00076 for (int cell=0; cell<result.dimension(0); cell++) { 00077 for (int pt=0; pt<result.dimension(1); pt++) { 00078 result(cell,pt) = pow(points(cell,pt,0), degree) - degree*(degree-1)*pow(points(cell,pt,0), degree-2); 00079 } 00080 } 00081 } 00082 } 00083 00085 void neumann(FieldContainer<double> & g_phi, const FieldContainer<double> & phi1, const FieldContainer<double> & phi2, int degree) { 00086 double g_at_one, g_at_minus_one; 00087 int num_fields; 00088 00089 if (degree == 0) { 00090 g_at_one = 0.0; 00091 g_at_minus_one = 0.0; 00092 } 00093 else { 00094 g_at_one = degree*pow(1.0, degree-1); 00095 g_at_minus_one = degree*pow(-1.0, degree-1); 00096 } 00097 00098 num_fields = phi1.dimension(0); 00099 00100 for (int i=0; i<num_fields; i++) { 00101 g_phi(0,i) = g_at_minus_one*phi1(i,0); 00102 g_phi(1,i) = g_at_one*phi2(i,0); 00103 } 00104 } 00105 00107 void u_exact(FieldContainer<double> & result, const FieldContainer<double> & points, int degree) { 00108 for (int cell=0; cell<result.dimension(0); cell++) { 00109 for (int pt=0; pt<result.dimension(1); pt++) { 00110 result(cell,pt) = pow(points(pt,0), degree); 00111 } 00112 } 00113 } 00114 00115 00116 00117 00118 int main(int argc, char *argv[]) { 00119 00120 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00121 00122 // This little trick lets us print to std::cout only if 00123 // a (dummy) command-line argument is provided. 00124 int iprint = argc - 1; 00125 Teuchos::RCP<std::ostream> outStream; 00126 Teuchos::oblackholestream bhs; // outputs nothing 00127 if (iprint > 0) 00128 outStream = Teuchos::rcp(&std::cout, false); 00129 else 00130 outStream = Teuchos::rcp(&bhs, false); 00131 00132 // Save the format state of the original std::cout. 00133 Teuchos::oblackholestream oldFormatState; 00134 oldFormatState.copyfmt(std::cout); 00135 00136 *outStream \ 00137 << "===============================================================================\n" \ 00138 << "| |\n" \ 00139 << "| Unit Test (Basis_HGRAD_LINE_Cn_FEM) |\n" \ 00140 << "| |\n" \ 00141 << "| 1) Patch test involving mass and stiffness matrices, |\n" \ 00142 << "| for the Neumann problem on a REFERENCE line: |\n" \ 00143 << "| |\n" \ 00144 << "| - u'' + u = f in (-1,1), u' = g at -1,1 |\n" \ 00145 << "| |\n" \ 00146 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00147 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00148 << "| Robert Kirby (robert.c.kirby@ttu.gov), |\n" \ 00149 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00150 << "| |\n" \ 00151 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00152 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00153 << "| |\n" \ 00154 << "===============================================================================\n"\ 00155 << "| TEST 1: Patch test |\n"\ 00156 << "===============================================================================\n"; 00157 00158 00159 int errorFlag = 0; 00160 double zero = 100*INTREPID_TOL; 00161 outStream -> precision(20); 00162 00163 00164 try { 00165 00166 int max_order = 10; // max total order of polynomial solution 00167 00168 // Define array containing points at which the solution is evaluated 00169 int numIntervals = 100; 00170 int numInterpPoints = numIntervals + 1; 00171 FieldContainer<double> interp_points(numInterpPoints, 1); 00172 for (int i=0; i<numInterpPoints; i++) { 00173 interp_points(i,0) = -1.0+(2.0*(double)i)/(double)numIntervals; 00174 } 00175 00176 DefaultCubatureFactory<double> cubFactory; // create factory 00177 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >()); // create cell topology 00178 00179 00180 for (int basis_order=1; basis_order <= max_order; basis_order++) { 00181 //create basis 00182 Teuchos::RCP<Basis<double,FieldContainer<double> > > lineBasis = 00183 Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<double,FieldContainer<double> >(basis_order,POINTTYPE_SPECTRAL) ); 00184 00185 for (int soln_order=1; soln_order <= basis_order; soln_order++) { 00186 // evaluate exact solution 00187 FieldContainer<double> exact_solution(1, numInterpPoints); 00188 u_exact(exact_solution, interp_points, soln_order); 00189 00190 int numFields = lineBasis->getCardinality(); 00191 00192 // create cubature 00193 Teuchos::RCP<Cubature<double> > lineCub = cubFactory.create(line, 2*basis_order-2); 00194 int numCubPoints = lineCub->getNumPoints(); 00195 int spaceDim = lineCub->getDimension(); 00196 00197 /* Computational arrays. */ 00198 FieldContainer<double> cub_points(numCubPoints, spaceDim); 00199 FieldContainer<double> cub_points_physical(1, numCubPoints, spaceDim); 00200 FieldContainer<double> cub_weights(numCubPoints); 00201 FieldContainer<double> cell_nodes(1, 2, spaceDim); 00202 FieldContainer<double> jacobian(1, numCubPoints, spaceDim, spaceDim); 00203 FieldContainer<double> jacobian_inv(1, numCubPoints, spaceDim, spaceDim); 00204 FieldContainer<double> jacobian_det(1, numCubPoints); 00205 FieldContainer<double> weighted_measure(1, numCubPoints); 00206 00207 FieldContainer<double> value_of_basis_at_cub_points(numFields, numCubPoints); 00208 FieldContainer<double> transformed_value_of_basis_at_cub_points(1, numFields, numCubPoints); 00209 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points(1, numFields, numCubPoints); 00210 FieldContainer<double> grad_of_basis_at_cub_points(numFields, numCubPoints, spaceDim); 00211 FieldContainer<double> transformed_grad_of_basis_at_cub_points(1, numFields, numCubPoints, spaceDim); 00212 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points(1, numFields, numCubPoints, spaceDim); 00213 FieldContainer<double> fe_matrix(1, numFields, numFields); 00214 00215 FieldContainer<double> rhs_at_cub_points_physical(1, numCubPoints); 00216 FieldContainer<double> rhs_and_soln_vector(1, numFields); 00217 00218 FieldContainer<double> one_point(1, 1); 00219 FieldContainer<double> value_of_basis_at_one(numFields, 1); 00220 FieldContainer<double> value_of_basis_at_minusone(numFields, 1); 00221 FieldContainer<double> bc_neumann(2, numFields); 00222 00223 FieldContainer<double> value_of_basis_at_interp_points(numFields, numInterpPoints); 00224 FieldContainer<double> transformed_value_of_basis_at_interp_points(1, numFields, numInterpPoints); 00225 FieldContainer<double> interpolant(1, numInterpPoints); 00226 00227 FieldContainer<int> ipiv(numFields); 00228 00229 /******************* START COMPUTATION ***********************/ 00230 00231 // get cubature points and weights 00232 lineCub->getCubature(cub_points, cub_weights); 00233 00234 // fill cell vertex array 00235 cell_nodes(0, 0, 0) = -1.0; 00236 cell_nodes(0, 1, 0) = 1.0; 00237 00238 // compute geometric cell information 00239 CellTools<double>::setJacobian(jacobian, cub_points, cell_nodes, line); 00240 CellTools<double>::setJacobianInv(jacobian_inv, jacobian); 00241 CellTools<double>::setJacobianDet(jacobian_det, jacobian); 00242 00243 // compute weighted measure 00244 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure, jacobian_det, cub_weights); 00245 00247 // Computing mass matrices: 00248 // tabulate values of basis functions at (reference) cubature points 00249 lineBasis->getValues(value_of_basis_at_cub_points, cub_points, OPERATOR_VALUE); 00250 00251 // transform values of basis functions 00252 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points, 00253 value_of_basis_at_cub_points); 00254 00255 // multiply with weighted measure 00256 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points, 00257 weighted_measure, 00258 transformed_value_of_basis_at_cub_points); 00259 00260 // compute mass matrices 00261 FunctionSpaceTools::integrate<double>(fe_matrix, 00262 transformed_value_of_basis_at_cub_points, 00263 weighted_transformed_value_of_basis_at_cub_points, 00264 COMP_CPP); 00266 00268 // Computing stiffness matrices: 00269 // tabulate gradients of basis functions at (reference) cubature points 00270 lineBasis->getValues(grad_of_basis_at_cub_points, cub_points, OPERATOR_GRAD); 00271 00272 // transform gradients of basis functions 00273 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points, 00274 jacobian_inv, 00275 grad_of_basis_at_cub_points); 00276 00277 // multiply with weighted measure 00278 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points, 00279 weighted_measure, 00280 transformed_grad_of_basis_at_cub_points); 00281 00282 // compute stiffness matrices and sum into fe_matrix 00283 FunctionSpaceTools::integrate<double>(fe_matrix, 00284 transformed_grad_of_basis_at_cub_points, 00285 weighted_transformed_grad_of_basis_at_cub_points, 00286 COMP_CPP, 00287 true); 00289 00291 // Computing RHS contributions: 00292 // map (reference) cubature points to physical space 00293 CellTools<double>::mapToPhysicalFrame(cub_points_physical, cub_points, cell_nodes, line); 00294 00295 // evaluate rhs function 00296 rhsFunc(rhs_at_cub_points_physical, cub_points_physical, soln_order); 00297 00298 // compute rhs 00299 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector, 00300 rhs_at_cub_points_physical, 00301 weighted_transformed_value_of_basis_at_cub_points, 00302 COMP_CPP); 00303 00304 // compute neumann b.c. contributions and adjust rhs 00305 one_point(0,0) = 1.0; lineBasis->getValues(value_of_basis_at_one, one_point, OPERATOR_VALUE); 00306 one_point(0,0) = -1.0; lineBasis->getValues(value_of_basis_at_minusone, one_point, OPERATOR_VALUE); 00307 neumann(bc_neumann, value_of_basis_at_minusone, value_of_basis_at_one, soln_order); 00308 for (int i=0; i<numFields; i++) { 00309 rhs_and_soln_vector(0, i) -= bc_neumann(0, i); 00310 rhs_and_soln_vector(0, i) += bc_neumann(1, i); 00311 } 00313 00315 // Solution of linear system: 00316 int info = 0; 00317 Teuchos::LAPACK<int, double> solver; 00318 //solver.GESV(numRows, 1, &fe_mat(0,0), numRows, &ipiv(0), &fe_vec(0), numRows, &info); 00319 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info); 00321 00323 // Building interpolant: 00324 // evaluate basis at interpolation points 00325 lineBasis->getValues(value_of_basis_at_interp_points, interp_points, OPERATOR_VALUE); 00326 // transform values of basis functions 00327 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points, 00328 value_of_basis_at_interp_points); 00329 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points); 00331 00332 /******************* END COMPUTATION ***********************/ 00333 00334 RealSpaceTools<double>::subtract(interpolant, exact_solution); 00335 00336 *outStream << "\nNorm-2 difference between exact solution polynomial of order " 00337 << soln_order << " and finite element interpolant of order " << basis_order << ": " 00338 << RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) << "\n"; 00339 00340 if (RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) > zero) { 00341 *outStream << "\n\nPatch test failed for solution polynomial order " 00342 << soln_order << " and basis order " << basis_order << "\n\n"; 00343 errorFlag++; 00344 } 00345 00346 } // end for basis_order 00347 00348 } // end for soln_order 00349 00350 } 00351 // Catch unexpected errors 00352 catch (std::logic_error err) { 00353 *outStream << err.what() << "\n\n"; 00354 errorFlag = -1000; 00355 }; 00356 00357 if (errorFlag != 0) 00358 std::cout << "End Result: TEST FAILED\n"; 00359 else 00360 std::cout << "End Result: TEST PASSED\n"; 00361 00362 // reset format state of std::cout 00363 std::cout.copyfmt(oldFormatState); 00364 00365 return errorFlag; 00366 }
1.7.4