00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00036 #include "Intrepid_FieldContainer.hpp"
00037 #include "Intrepid_HGRAD_TRI_C2_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 "Teuchos_oblackholestream.hpp"
00044 #include "Teuchos_RCP.hpp"
00045 #include "Teuchos_GlobalMPISession.hpp"
00046 #include "Teuchos_SerialDenseMatrix.hpp"
00047 #include "Teuchos_SerialDenseVector.hpp"
00048 #include "Teuchos_LAPACK.hpp"
00049
00050 using namespace std;
00051 using namespace Intrepid;
00052
00053 void rhsFunc(FieldContainer<double> &, const FieldContainer<double> &, int, int);
00054 void neumann(FieldContainer<double> & ,
00055 const FieldContainer<double> & ,
00056 const FieldContainer<double> & ,
00057 const shards::CellTopology & ,
00058 int, int, int);
00059 void u_exact(FieldContainer<double> &, const FieldContainer<double> &, int, int);
00060
00062 void rhsFunc(FieldContainer<double> & result,
00063 const FieldContainer<double> & points,
00064 int xd,
00065 int yd) {
00066
00067 int x = 0, y = 1;
00068
00069
00070 if (xd > 1) {
00071 for (int cell=0; cell<result.dimension(0); cell++) {
00072 for (int pt=0; pt<result.dimension(1); pt++) {
00073 result(cell,pt) = - xd*(xd-1)*std::pow(points(cell,pt,x), xd-2) * std::pow(points(cell,pt,y), yd);
00074 }
00075 }
00076 }
00077
00078
00079 if (yd > 1) {
00080 for (int cell=0; cell<result.dimension(0); cell++) {
00081 for (int pt=0; pt<result.dimension(1); pt++) {
00082 result(cell,pt) -= yd*(yd-1)*std::pow(points(cell,pt,y), yd-2) * std::pow(points(cell,pt,x), xd);
00083 }
00084 }
00085 }
00086
00087
00088 for (int cell=0; cell<result.dimension(0); cell++) {
00089 for (int pt=0; pt<result.dimension(1); pt++) {
00090 result(cell,pt) += std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
00091 }
00092 }
00093
00094 }
00095
00096
00098 void neumann(FieldContainer<double> & result,
00099 const FieldContainer<double> & points,
00100 const FieldContainer<double> & jacs,
00101 const shards::CellTopology & parentCell,
00102 int sideOrdinal, int xd, int yd) {
00103
00104 int x = 0, y = 1;
00105
00106 int numCells = result.dimension(0);
00107 int numPoints = result.dimension(1);
00108
00109 FieldContainer<double> grad_u(numCells, numPoints, 2);
00110 FieldContainer<double> side_normals(numCells, numPoints, 2);
00111 FieldContainer<double> normal_lengths(numCells, numPoints);
00112
00113
00114 if (xd > 0) {
00115 for (int cell=0; cell<numCells; cell++) {
00116 for (int pt=0; pt<numPoints; pt++) {
00117 grad_u(cell,pt,x) = xd*std::pow(points(cell,pt,x), xd-1) * std::pow(points(cell,pt,y), yd);
00118 }
00119 }
00120 }
00121
00122
00123 if (yd > 0) {
00124 for (int cell=0; cell<numCells; cell++) {
00125 for (int pt=0; pt<numPoints; pt++) {
00126 grad_u(cell,pt,y) = yd*std::pow(points(cell,pt,y), yd-1) * std::pow(points(cell,pt,x), xd);
00127 }
00128 }
00129 }
00130
00131 CellTools<double>::getPhysicalSideNormals(side_normals, jacs, sideOrdinal, parentCell);
00132
00133
00134 RealSpaceTools<double>::vectorNorm(normal_lengths, side_normals, NORM_TWO);
00135 FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals, true);
00136
00137 FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals);
00138
00139 }
00140
00142 void u_exact(FieldContainer<double> & result, const FieldContainer<double> & points, int xd, int yd) {
00143 int x = 0, y = 1;
00144 for (int cell=0; cell<result.dimension(0); cell++) {
00145 for (int pt=0; pt<result.dimension(1); pt++) {
00146 result(cell,pt) = std::pow(points(pt,x), xd)*std::pow(points(pt,y), yd);
00147 }
00148 }
00149 }
00150
00151
00152
00153
00154 int main(int argc, char *argv[]) {
00155
00156 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00157
00158
00159
00160 int iprint = argc - 1;
00161 Teuchos::RCP<std::ostream> outStream;
00162 Teuchos::oblackholestream bhs;
00163 if (iprint > 0)
00164 outStream = Teuchos::rcp(&std::cout, false);
00165 else
00166 outStream = Teuchos::rcp(&bhs, false);
00167
00168
00169 Teuchos::oblackholestream oldFormatState;
00170 oldFormatState.copyfmt(std::cout);
00171
00172 *outStream \
00173 << "===============================================================================\n" \
00174 << "| |\n" \
00175 << "| Unit Test (Basis_HGRAD_TRI_C2_FEM) |\n" \
00176 << "| |\n" \
00177 << "| 1) Patch test involving mass and stiffness matrices, |\n" \
00178 << "| for the Neumann problem on a triangular patch |\n" \
00179 << "| Omega with boundary Gamma. |\n" \
00180 << "| |\n" \
00181 << "| - div (grad u) + u = f in Omega, (grad u) . n = g on Gamma |\n" \
00182 << "| |\n" \
00183 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
00184 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
00185 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
00186 << "| |\n" \
00187 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
00188 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
00189 << "| |\n" \
00190 << "===============================================================================\n"\
00191 << "| TEST 1: Patch test |\n"\
00192 << "===============================================================================\n";
00193
00194
00195 int errorFlag = 0;
00196
00197 outStream -> precision(16);
00198
00199
00200 try {
00201
00202 int max_order = 2;
00203 DefaultCubatureFactory<double> cubFactory;
00204 shards::CellTopology cell(shards::getCellTopologyData< shards::Triangle<> >());
00205 shards::CellTopology side(shards::getCellTopologyData< shards::Line<> >());
00206 int cellDim = cell.getDimension();
00207 int sideDim = side.getDimension();
00208
00209
00210 int numIntervals = 10;
00211 int numInterpPoints = ((numIntervals + 1)*(numIntervals + 2))/2;
00212 FieldContainer<double> interp_points_ref(numInterpPoints, 2);
00213 int counter = 0;
00214 for (int j=0; j<=numIntervals; j++) {
00215 for (int i=0; i<=numIntervals; i++) {
00216 if (i <= numIntervals-j) {
00217 interp_points_ref(counter,0) = i*(1.0/numIntervals);
00218 interp_points_ref(counter,1) = j*(1.0/numIntervals);
00219 counter++;
00220 }
00221 }
00222 }
00223
00224
00225 FieldContainer<double> cell_nodes(1, 3, cellDim);
00226
00227 cell_nodes(0, 0, 0) = 0.1;
00228 cell_nodes(0, 0, 1) = -0.1;
00229 cell_nodes(0, 1, 0) = 1.1;
00230 cell_nodes(0, 1, 1) = -0.1;
00231 cell_nodes(0, 2, 0) = 0.1;
00232 cell_nodes(0, 2, 1) = 0.9;
00233
00234
00235
00236
00237
00238
00239
00240
00241 FieldContainer<double> interp_points(1, numInterpPoints, cellDim);
00242 CellTools<double>::mapToPhysicalFrame(interp_points, interp_points_ref, cell_nodes, cell);
00243 interp_points.resize(numInterpPoints, cellDim);
00244
00245 for (int x_order=0; x_order <= max_order; x_order++) {
00246 for (int y_order=0; y_order <= max_order-x_order; y_order++) {
00247
00248
00249 FieldContainer<double> exact_solution(1, numInterpPoints);
00250 u_exact(exact_solution, interp_points, x_order, y_order);
00251
00252 int basis_order = 2;
00253
00254
00255 double zero = basis_order*basis_order*100*INTREPID_TOL;
00256
00257
00258 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
00259 Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM<double,FieldContainer<double> >() );
00260 int numFields = basis->getCardinality();
00261
00262
00263 Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*basis_order);
00264 Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*basis_order);
00265 int numCubPointsCell = cellCub->getNumPoints();
00266 int numCubPointsSide = sideCub->getNumPoints();
00267
00268
00269
00270 FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim);
00271 FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim);
00272 FieldContainer<double> cub_weights_cell(numCubPointsCell);
00273 FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim);
00274 FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim);
00275 FieldContainer<double> jacobian_det_cell(1, numCubPointsCell);
00276 FieldContainer<double> weighted_measure_cell(1, numCubPointsCell);
00277
00278 FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell);
00279 FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00280 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00281 FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim);
00282 FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00283 FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00284 FieldContainer<double> fe_matrix(1, numFields, numFields);
00285
00286 FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell);
00287 FieldContainer<double> rhs_and_soln_vector(1, numFields);
00288
00289
00290 unsigned numSides = 3;
00291 FieldContainer<double> cub_points_side(numCubPointsSide, sideDim);
00292 FieldContainer<double> cub_weights_side(numCubPointsSide);
00293 FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim);
00294 FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim);
00295 FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim);
00296 FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide);
00297 FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide);
00298
00299 FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide);
00300 FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00301 FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00302 FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide);
00303 FieldContainer<double> neumann_fields_per_side(1, numFields);
00304
00305
00306 FieldContainer<double> value_of_basis_at_interp_points(numFields, numInterpPoints);
00307 FieldContainer<double> transformed_value_of_basis_at_interp_points(1, numFields, numInterpPoints);
00308 FieldContainer<double> interpolant(1, numInterpPoints);
00309
00310 FieldContainer<int> ipiv(numFields);
00311
00312
00313
00314
00315
00316
00317 cellCub->getCubature(cub_points_cell, cub_weights_cell);
00318
00319
00320 CellTools<double>::setJacobian(jacobian_cell, cub_points_cell, cell_nodes, cell);
00321 CellTools<double>::setJacobianInv(jacobian_inv_cell, jacobian_cell);
00322 CellTools<double>::setJacobianDet(jacobian_det_cell, jacobian_cell);
00323
00324
00325 FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
00326
00328
00329
00330 basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
00331
00332
00333 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
00334 value_of_basis_at_cub_points_cell);
00335
00336
00337 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
00338 weighted_measure_cell,
00339 transformed_value_of_basis_at_cub_points_cell);
00340
00341
00342 FunctionSpaceTools::integrate<double>(fe_matrix,
00343 transformed_value_of_basis_at_cub_points_cell,
00344 weighted_transformed_value_of_basis_at_cub_points_cell,
00345 COMP_BLAS);
00347
00349
00350
00351 basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
00352
00353
00354 FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
00355 jacobian_inv_cell,
00356 grad_of_basis_at_cub_points_cell);
00357
00358
00359 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
00360 weighted_measure_cell,
00361 transformed_grad_of_basis_at_cub_points_cell);
00362
00363
00364 FunctionSpaceTools::integrate<double>(fe_matrix,
00365 transformed_grad_of_basis_at_cub_points_cell,
00366 weighted_transformed_grad_of_basis_at_cub_points_cell,
00367 COMP_BLAS,
00368 true);
00370
00372
00373
00374 CellTools<double>::mapToPhysicalFrame(cub_points_cell_physical, cub_points_cell, cell_nodes, cell);
00375
00376
00377 rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order);
00378
00379
00380 FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
00381 rhs_at_cub_points_cell_physical,
00382 weighted_transformed_value_of_basis_at_cub_points_cell,
00383 COMP_BLAS);
00384
00385
00386 sideCub->getCubature(cub_points_side, cub_weights_side);
00387 for (unsigned i=0; i<numSides; i++) {
00388
00389 CellTools<double>::mapToReferenceSubcell(cub_points_side_refcell, cub_points_side, sideDim, (int)i, cell);
00390 CellTools<double>::setJacobian(jacobian_side_refcell, cub_points_side_refcell, cell_nodes, cell);
00391 CellTools<double>::setJacobianDet(jacobian_det_side_refcell, jacobian_side_refcell);
00392
00393
00394 FunctionSpaceTools::computeEdgeMeasure<double>(weighted_measure_side_refcell,
00395 jacobian_side_refcell,
00396 cub_weights_side,
00397 i,
00398 cell);
00399
00400
00401 basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
00402
00403 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
00404 value_of_basis_at_cub_points_side_refcell);
00405
00406
00407 FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00408 weighted_measure_side_refcell,
00409 transformed_value_of_basis_at_cub_points_side_refcell);
00410
00411
00412
00413 CellTools<double>::mapToPhysicalFrame(cub_points_side_physical, cub_points_side_refcell, cell_nodes, cell);
00414
00415 neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
00416 cell, (int)i, x_order, y_order);
00417
00418 FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
00419 neumann_data_at_cub_points_side_physical,
00420 weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00421 COMP_BLAS);
00422
00423
00424 RealSpaceTools<double>::add(rhs_and_soln_vector, neumann_fields_per_side);;
00425 }
00427
00429
00430 int info = 0;
00431 Teuchos::LAPACK<int, double> solver;
00432 solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
00434
00436
00437
00438 basis->getValues(value_of_basis_at_interp_points, interp_points_ref, OPERATOR_VALUE);
00439
00440 FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points,
00441 value_of_basis_at_interp_points);
00442 FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points);
00444
00445
00446
00447 RealSpaceTools<double>::subtract(interpolant, exact_solution);
00448
00449 *outStream << "\nRelative norm-2 error between exact solution polynomial of order ("
00450 << x_order << ", " << y_order << ") and finite element interpolant of order " << basis_order << ": "
00451 << RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00452 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) << "\n";
00453
00454 if (RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00455 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) > zero) {
00456 *outStream << "\n\nPatch test failed for solution polynomial order ("
00457 << x_order << ", " << y_order << ") and basis order " << basis_order << "\n\n";
00458 errorFlag++;
00459 }
00460 }
00461 }
00462
00463 }
00464
00465 catch (std::logic_error err) {
00466 *outStream << err.what() << "\n\n";
00467 errorFlag = -1000;
00468 };
00469
00470 if (errorFlag != 0)
00471 std::cout << "End Result: TEST FAILED\n";
00472 else
00473 std::cout << "End Result: TEST PASSED\n";
00474
00475
00476 std::cout.copyfmt(oldFormatState);
00477
00478 return errorFlag;
00479 }