Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Basis/HGRAD_TRI_Cn_FEM/test_02.cpp
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_TRI_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 "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   // second x-derivatives of u
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   // second y-derivatives of u
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   // add u
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   // first x-derivatives of u
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   // first y-derivatives of u
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   // scale normals
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   // This little trick lets us print to std::cout only if
00159   // a (dummy) command-line argument is provided.
00160   int iprint     = argc - 1;
00161   Teuchos::RCP<std::ostream> outStream;
00162   Teuchos::oblackholestream bhs; // outputs nothing
00163   if (iprint > 0)
00164     outStream = Teuchos::rcp(&std::cout, false);
00165   else
00166     outStream = Teuchos::rcp(&bhs, false);
00167 
00168   // Save the format state of the original std::cout.
00169   Teuchos::oblackholestream oldFormatState;
00170   oldFormatState.copyfmt(std::cout);
00171 
00172   *outStream \
00173     << "===============================================================================\n" \
00174     << "|                                                                             |\n" \
00175     << "|                    Unit Test (Basis_HGRAD_TRI_Cn_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 = 10;                                                              // max total order of polynomial solution
00203     DefaultCubatureFactory<double>  cubFactory;                                      // create cubature factory
00204     shards::CellTopology cell(shards::getCellTopologyData< shards::Triangle<> >());  // create parent cell topology
00205     shards::CellTopology side(shards::getCellTopologyData< shards::Line<> >());      // create relevant subcell (side) topology
00206     int cellDim = cell.getDimension();
00207     int sideDim = side.getDimension();
00208 
00209     // Define array containing points at which the solution is evaluated, in reference cell.
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     /* Parent cell definition. */
00225     FieldContainer<double> cell_nodes(1, 3, cellDim);
00226     // Perturbed reference triangle.
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     // Reference triangle. 
00234     /*cell_nodes(0, 0, 0) = 0.0;
00235     cell_nodes(0, 0, 1) = 0.0;
00236     cell_nodes(0, 1, 0) = 1.0;
00237     cell_nodes(0, 1, 1) = 0.0;
00238     cell_nodes(0, 2, 0) = 0.0;
00239     cell_nodes(0, 2, 1) = 1.0;*/
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     // we'll test two types of bases
00246     EPointType pointtype[] = {POINTTYPE_EQUISPACED, POINTTYPE_WARPBLEND};
00247     for (int ptype=0; ptype < 2; ptype++) {
00248 
00249       *outStream << "\nTesting bases with " << EPointTypeToString(pointtype[ptype]) << ":\n";
00250 
00251       for (int x_order=0; x_order <= max_order; x_order++) {
00252         for (int y_order=0; y_order <= max_order-x_order; y_order++) {
00253 
00254           // evaluate exact solution
00255           FieldContainer<double> exact_solution(1, numInterpPoints);
00256           u_exact(exact_solution, interp_points, x_order, y_order);
00257 
00258           int total_order = std::max(x_order + y_order, 1);
00259 
00260           for (int basis_order=total_order; basis_order <= max_order; basis_order++) {
00261 
00262             // set test tolerance
00263             double zero = basis_order*basis_order*100*INTREPID_TOL;
00264 
00265             //create basis
00266             Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
00267               Teuchos::rcp(new Basis_HGRAD_TRI_Cn_FEM<double,FieldContainer<double> >(basis_order, pointtype[ptype]) );
00268             int numFields = basis->getCardinality();
00269 
00270             // create cubatures
00271             Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*basis_order);
00272             Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*basis_order);
00273             int numCubPointsCell = cellCub->getNumPoints();
00274             int numCubPointsSide = sideCub->getNumPoints();
00275 
00276             /* Computational arrays. */
00277             /* Section 1: Related to parent cell integration. */
00278             FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim);
00279             FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim);
00280             FieldContainer<double> cub_weights_cell(numCubPointsCell);
00281             FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim);
00282             FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim);
00283             FieldContainer<double> jacobian_det_cell(1, numCubPointsCell);
00284             FieldContainer<double> weighted_measure_cell(1, numCubPointsCell);
00285 
00286             FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell);
00287             FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00288             FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00289             FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim);
00290             FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00291             FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00292             FieldContainer<double> fe_matrix(1, numFields, numFields);
00293 
00294             FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell);
00295             FieldContainer<double> rhs_and_soln_vector(1, numFields);
00296 
00297             /* Section 2: Related to subcell (side) integration. */
00298             unsigned numSides = 3;
00299             FieldContainer<double> cub_points_side(numCubPointsSide, sideDim);
00300             FieldContainer<double> cub_weights_side(numCubPointsSide);
00301             FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim);
00302             FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim);
00303             FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim);
00304             FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide);
00305             FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide);
00306 
00307             FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide);
00308             FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00309             FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00310             FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide);
00311             FieldContainer<double> neumann_fields_per_side(1, numFields);
00312 
00313             /* Section 3: Related to global interpolant. */
00314             FieldContainer<double> value_of_basis_at_interp_points(numFields, numInterpPoints);
00315             FieldContainer<double> transformed_value_of_basis_at_interp_points(1, numFields, numInterpPoints);
00316             FieldContainer<double> interpolant(1, numInterpPoints);
00317 
00318             FieldContainer<int> ipiv(numFields);
00319 
00320 
00321 
00322             /******************* START COMPUTATION ***********************/
00323 
00324             // get cubature points and weights
00325             cellCub->getCubature(cub_points_cell, cub_weights_cell);
00326 
00327             // compute geometric cell information
00328             CellTools<double>::setJacobian(jacobian_cell, cub_points_cell, cell_nodes, cell);
00329             CellTools<double>::setJacobianInv(jacobian_inv_cell, jacobian_cell);
00330             CellTools<double>::setJacobianDet(jacobian_det_cell, jacobian_cell);
00331 
00332             // compute weighted measure
00333             FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
00334 
00336             // Computing mass matrices:
00337             // tabulate values of basis functions at (reference) cubature points
00338             basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
00339 
00340             // transform values of basis functions
00341             FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
00342                                                             value_of_basis_at_cub_points_cell);
00343 
00344             // multiply with weighted measure
00345             FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
00346                                                         weighted_measure_cell,
00347                                                         transformed_value_of_basis_at_cub_points_cell);
00348 
00349             // compute mass matrices
00350             FunctionSpaceTools::integrate<double>(fe_matrix,
00351                                                   transformed_value_of_basis_at_cub_points_cell,
00352                                                   weighted_transformed_value_of_basis_at_cub_points_cell,
00353                                                   COMP_BLAS);
00355 
00357             // Computing stiffness matrices:
00358             // tabulate gradients of basis functions at (reference) cubature points
00359             basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
00360 
00361             // transform gradients of basis functions
00362             FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
00363                                                            jacobian_inv_cell,
00364                                                            grad_of_basis_at_cub_points_cell);
00365 
00366             // multiply with weighted measure
00367             FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
00368                                                         weighted_measure_cell,
00369                                                         transformed_grad_of_basis_at_cub_points_cell);
00370 
00371             // compute stiffness matrices and sum into fe_matrix
00372             FunctionSpaceTools::integrate<double>(fe_matrix,
00373                                                   transformed_grad_of_basis_at_cub_points_cell,
00374                                                   weighted_transformed_grad_of_basis_at_cub_points_cell,
00375                                                   COMP_BLAS,
00376                                                   true);
00378 
00380             // Computing RHS contributions:
00381             // map cell (reference) cubature points to physical space
00382             CellTools<double>::mapToPhysicalFrame(cub_points_cell_physical, cub_points_cell, cell_nodes, cell);
00383 
00384             // evaluate rhs function
00385             rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order);
00386 
00387             // compute rhs
00388             FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
00389                                                   rhs_at_cub_points_cell_physical,
00390                                                   weighted_transformed_value_of_basis_at_cub_points_cell,
00391                                                   COMP_BLAS);
00392 
00393             // compute neumann b.c. contributions and adjust rhs
00394             sideCub->getCubature(cub_points_side, cub_weights_side);
00395             for (unsigned i=0; i<numSides; i++) {
00396               // compute geometric cell information
00397               CellTools<double>::mapToReferenceSubcell(cub_points_side_refcell, cub_points_side, sideDim, (int)i, cell);
00398               CellTools<double>::setJacobian(jacobian_side_refcell, cub_points_side_refcell, cell_nodes, cell);
00399               CellTools<double>::setJacobianDet(jacobian_det_side_refcell, jacobian_side_refcell);
00400 
00401               // compute weighted edge measure
00402               FunctionSpaceTools::computeEdgeMeasure<double>(weighted_measure_side_refcell,
00403                                                              jacobian_side_refcell,
00404                                                              cub_weights_side,
00405                                                              i,
00406                                                              cell);
00407 
00408               // tabulate values of basis functions at side cubature points, in the reference parent cell domain
00409               basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
00410               // transform
00411               FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
00412                                                               value_of_basis_at_cub_points_side_refcell);
00413 
00414               // multiply with weighted measure
00415               FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00416                                                           weighted_measure_side_refcell,
00417                                                           transformed_value_of_basis_at_cub_points_side_refcell);
00418 
00419               // compute Neumann data
00420               // map side cubature points in reference parent cell domain to physical space
00421               CellTools<double>::mapToPhysicalFrame(cub_points_side_physical, cub_points_side_refcell, cell_nodes, cell);
00422               // now compute data
00423               neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
00424                       cell, (int)i, x_order, y_order);
00425 
00426               FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
00427                                                     neumann_data_at_cub_points_side_physical,
00428                                                     weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00429                                                     COMP_BLAS);
00430 
00431               // adjust RHS
00432               RealSpaceTools<double>::add(rhs_and_soln_vector, neumann_fields_per_side);;
00433             }
00435 
00437             // Solution of linear system:
00438             int info = 0;
00439             Teuchos::LAPACK<int, double> solver;
00440             solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
00442 
00444             // Building interpolant:
00445             // evaluate basis at interpolation points
00446             basis->getValues(value_of_basis_at_interp_points, interp_points_ref, OPERATOR_VALUE);
00447             // transform values of basis functions
00448             FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points,
00449                                                             value_of_basis_at_interp_points);
00450             FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points);
00452 
00453             /******************* END COMPUTATION ***********************/
00454         
00455             RealSpaceTools<double>::subtract(interpolant, exact_solution);
00456 
00457             *outStream << "\nRelative norm-2 error between exact solution polynomial of order ("
00458                        << x_order << ", " << y_order << ") and finite element interpolant of order " << basis_order << ": "
00459                        << RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00460                           RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) << "\n";
00461 
00462             if (RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00463                 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) > zero) {
00464               *outStream << "\n\nPatch test failed for solution polynomial order ("
00465                          << x_order << ", " << y_order << ") and basis order " << basis_order << "\n\n";
00466               errorFlag++;
00467             }
00468           } // end for basis_order
00469         } // end for y_order
00470       } // end for x_order
00471     } // end for ptype
00472 
00473   }
00474   // Catch unexpected errors
00475   catch (std::logic_error err) {
00476     *outStream << err.what() << "\n\n";
00477     errorFlag = -1000;
00478   };
00479 
00480   if (errorFlag != 0)
00481     std::cout << "End Result: TEST FAILED\n";
00482   else
00483     std::cout << "End Result: TEST PASSED\n";
00484 
00485   // reset format state of std::cout
00486   std::cout.copyfmt(oldFormatState);
00487 
00488   return errorFlag;
00489 }