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