Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/test/Discretization/Basis/HGRAD_HEX_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_HEX_Cn_FEM.hpp"
00038 #include "Intrepid_PointTools.hpp"
00039 #include "Intrepid_DefaultCubatureFactory.hpp"
00040 #include "Intrepid_RealSpaceTools.hpp"
00041 #include "Intrepid_ArrayTools.hpp"
00042 #include "Intrepid_FunctionSpaceTools.hpp"
00043 #include "Intrepid_CellTools.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, int, int);
00055 void neumann(FieldContainer<double>       & ,
00056              const FieldContainer<double> & ,
00057              const FieldContainer<double> & ,
00058              const shards::CellTopology   & ,
00059              int, int, int, int);
00060 void u_exact(FieldContainer<double> &, const FieldContainer<double> &, int, int, int);
00061 
00063 void rhsFunc(FieldContainer<double> & result,
00064              const FieldContainer<double> & points,
00065              int xd,
00066              int yd,
00067              int zd) {
00068 
00069   int x = 0, y = 1, z = 2;
00070 
00071   // second x-derivatives of u
00072   if (xd > 1) {
00073     for (int cell=0; cell<result.dimension(0); cell++) {
00074       for (int pt=0; pt<result.dimension(1); pt++) {
00075         result(cell,pt) = - xd*(xd-1)*std::pow(points(cell,pt,x), xd-2) *
00076                             std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
00077       }
00078     }
00079   }
00080 
00081   // second y-derivatives of u
00082   if (yd > 1) {
00083     for (int cell=0; cell<result.dimension(0); cell++) {
00084       for (int pt=0; pt<result.dimension(1); pt++) {
00085         result(cell,pt) -=  yd*(yd-1)*std::pow(points(cell,pt,y), yd-2) *
00086                             std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd);
00087       }
00088     }
00089   }
00090 
00091   // second z-derivatives of u
00092   if (zd > 1) {
00093     for (int cell=0; cell<result.dimension(0); cell++) {
00094       for (int pt=0; pt<result.dimension(1); pt++) {
00095         result(cell,pt) -=  zd*(zd-1)*std::pow(points(cell,pt,z), zd-2) *
00096                             std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
00097       }
00098     }
00099   }
00100 
00101   // add u
00102   for (int cell=0; cell<result.dimension(0); cell++) {
00103     for (int pt=0; pt<result.dimension(1); pt++) {
00104       result(cell,pt) +=  std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
00105     }
00106   }
00107 
00108 }
00109 
00110 
00112 void neumann(FieldContainer<double>       & result,
00113              const FieldContainer<double> & points,
00114              const FieldContainer<double> & jacs,
00115              const shards::CellTopology   & parentCell,
00116              int sideOrdinal, int xd, int yd, int zd) {
00117 
00118   int x = 0, y = 1, z = 2;
00119 
00120   int numCells  = result.dimension(0);
00121   int numPoints = result.dimension(1);
00122 
00123   FieldContainer<double> grad_u(numCells, numPoints, 3);
00124   FieldContainer<double> side_normals(numCells, numPoints, 3);
00125   FieldContainer<double> normal_lengths(numCells, numPoints);
00126 
00127   // first x-derivatives of u
00128   if (xd > 0) {
00129     for (int cell=0; cell<numCells; cell++) {
00130       for (int pt=0; pt<numPoints; pt++) {
00131         grad_u(cell,pt,x) = xd*std::pow(points(cell,pt,x), xd-1) *
00132                             std::pow(points(cell,pt,y), yd) * std::pow(points(cell,pt,z), zd);
00133       }
00134     }
00135   }
00136 
00137   // first y-derivatives of u
00138   if (yd > 0) {
00139     for (int cell=0; cell<numCells; cell++) {
00140       for (int pt=0; pt<numPoints; pt++) {
00141         grad_u(cell,pt,y) = yd*std::pow(points(cell,pt,y), yd-1) *
00142                             std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,z), zd);
00143       }
00144     }
00145   }
00146 
00147   // first z-derivatives of u
00148   if (zd > 0) {
00149     for (int cell=0; cell<numCells; cell++) {
00150       for (int pt=0; pt<numPoints; pt++) {
00151         grad_u(cell,pt,z) = zd*std::pow(points(cell,pt,z), zd-1) *
00152                             std::pow(points(cell,pt,x), xd) * std::pow(points(cell,pt,y), yd);
00153       }
00154     }
00155   }
00156   
00157   CellTools<double>::getPhysicalSideNormals(side_normals, jacs, sideOrdinal, parentCell);
00158 
00159   // scale normals
00160   RealSpaceTools<double>::vectorNorm(normal_lengths, side_normals, NORM_TWO);
00161   FunctionSpaceTools::scalarMultiplyDataData<double>(side_normals, normal_lengths, side_normals, true); 
00162 
00163   FunctionSpaceTools::dotMultiplyDataData<double>(result, grad_u, side_normals);
00164 
00165 }
00166 
00168 void u_exact(FieldContainer<double> & result, const FieldContainer<double> & points, int xd, int yd, int zd) {
00169   int x = 0, y = 1, z = 2;
00170   for (int cell=0; cell<result.dimension(0); cell++) {
00171     for (int pt=0; pt<result.dimension(1); pt++) {
00172       result(cell,pt) = std::pow(points(pt,x), xd)*std::pow(points(pt,y), yd)*std::pow(points(pt,z), zd);
00173     }
00174   }
00175 }
00176 
00177 
00178 
00179 
00180 int main(int argc, char *argv[]) {
00181 
00182   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00183 
00184   // This little trick lets us print to std::cout only if
00185   // a (dummy) command-line argument is provided.
00186   int iprint     = argc - 1;
00187   Teuchos::RCP<std::ostream> outStream;
00188   Teuchos::oblackholestream bhs; // outputs nothing
00189   if (iprint > 0)
00190     outStream = Teuchos::rcp(&std::cout, false);
00191   else
00192     outStream = Teuchos::rcp(&bhs, false);
00193 
00194   // Save the format state of the original std::cout.
00195   Teuchos::oblackholestream oldFormatState;
00196   oldFormatState.copyfmt(std::cout);
00197 
00198   *outStream \
00199     << "===============================================================================\n" \
00200     << "|                                                                             |\n" \
00201     << "|                    Unit Test (Basis_HGRAD_HEX_Cn_FEM)                       |\n" \
00202     << "|                                                                             |\n" \
00203     << "|     1) Patch test involving mass and stiffness matrices,                    |\n" \
00204     << "|        for the Neumann problem on a physical parallelepiped                 |\n" \
00205     << "|        AND a reference hex Omega with boundary Gamma.                       |\n" \
00206     << "|                                                                             |\n" \
00207     << "|        - div (grad u) + u = f  in Omega,  (grad u) . n = g  on Gamma        |\n" \
00208     << "|                                                                             |\n" \
00209     << "|        For a generic parallelepiped, the basis recovers a complete          |\n" \
00210     << "|        polynomial space of order 2. On a (scaled and/or translated)         |\n" \
00211     << "|        reference hex, the basis recovers a complete tensor product          |\n" \
00212     << "|        space of order 1 (i.e. incl. cross terms, e.g. x^2*y^2*z^2).         |\n" \
00213     << "|                                                                             |\n" \
00214     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00215     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00216     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00217     << "|                                                                             |\n" \
00218     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00219     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00220     << "|                                                                             |\n" \
00221     << "===============================================================================\n"\
00222     << "| TEST 1: Patch test                                                          |\n"\
00223     << "===============================================================================\n";
00224 
00225   
00226   int errorFlag = 0;
00227 
00228   outStream -> precision(16);
00229 
00230 
00231   try {
00232 
00233     int max_order = 3;                                                                    // max total order of polynomial solution
00234     DefaultCubatureFactory<double>  cubFactory;                                           // create factory
00235     shards::CellTopology cell(shards::getCellTopologyData< shards::Hexahedron<> >());     // create parent cell topology
00236     shards::CellTopology side(shards::getCellTopologyData< shards::Quadrilateral<> >());  // create relevant subcell (side) topology
00237     shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());  // create relevant subcell (side) topology
00238     int cellDim = cell.getDimension();
00239     int sideDim = side.getDimension();
00240     unsigned numSides = 6;
00241 
00242     // Define array containing points at which the solution is evaluated, on the reference tet.
00243     int numIntervals = 10;
00244     int numInterpPoints = (numIntervals + 1)*(numIntervals + 1)*(numIntervals + 1);
00245     FieldContainer<double> interp_points_ref(numInterpPoints, 3);
00246     int counter = 0;
00247     for (int k=0; k<=numIntervals; k++) {
00248       for (int j=0; j<=numIntervals; j++) {
00249         for (int i=0; i<=numIntervals; i++) {
00250           interp_points_ref(counter,0) = i*(1.0/numIntervals)-1.0;
00251           interp_points_ref(counter,1) = j*(1.0/numIntervals)-1.0;
00252           interp_points_ref(counter,2) = k*(1.0/numIntervals)-1.0;
00253           counter++;
00254         }
00255       }
00256     }
00257 
00258     /* Parent cell definition. */
00259     FieldContainer<double> cell_nodes[2];
00260     cell_nodes[0].resize(1, 8, cellDim);
00261     cell_nodes[1].resize(1, 8, cellDim);
00262 
00263     // Generic parallelepiped.
00264     cell_nodes[0](0, 0, 0) = -5.0;
00265     cell_nodes[0](0, 0, 1) = -1.0;
00266     cell_nodes[0](0, 0, 2) = 0.0;
00267     cell_nodes[0](0, 1, 0) = 4.0;
00268     cell_nodes[0](0, 1, 1) = 1.0;
00269     cell_nodes[0](0, 1, 2) = 1.0;
00270     cell_nodes[0](0, 2, 0) = 8.0;
00271     cell_nodes[0](0, 2, 1) = 3.0;
00272     cell_nodes[0](0, 2, 2) = 1.0;
00273     cell_nodes[0](0, 3, 0) = -1.0;
00274     cell_nodes[0](0, 3, 1) = 1.0;
00275     cell_nodes[0](0, 3, 2) = 0.0;
00276     cell_nodes[0](0, 4, 0) = 5.0;
00277     cell_nodes[0](0, 4, 1) = 9.0;
00278     cell_nodes[0](0, 4, 2) = 1.0;
00279     cell_nodes[0](0, 5, 0) = 14.0;
00280     cell_nodes[0](0, 5, 1) = 11.0;
00281     cell_nodes[0](0, 5, 2) = 2.0;
00282     cell_nodes[0](0, 6, 0) = 18.0;
00283     cell_nodes[0](0, 6, 1) = 13.0;
00284     cell_nodes[0](0, 6, 2) = 2.0;
00285     cell_nodes[0](0, 7, 0) = 9.0;
00286     cell_nodes[0](0, 7, 1) = 11.0;
00287     cell_nodes[0](0, 7, 2) = 1.0;
00288     // Reference hex.
00289     cell_nodes[1](0, 0, 0) = -1.0;
00290     cell_nodes[1](0, 0, 1) = -1.0;
00291     cell_nodes[1](0, 0, 2) = -1.0;
00292     cell_nodes[1](0, 1, 0) = 1.0;
00293     cell_nodes[1](0, 1, 1) = -1.0;
00294     cell_nodes[1](0, 1, 2) = -1.0;
00295     cell_nodes[1](0, 2, 0) = 1.0;
00296     cell_nodes[1](0, 2, 1) = 1.0;
00297     cell_nodes[1](0, 2, 2) = -1.0;
00298     cell_nodes[1](0, 3, 0) = -1.0;
00299     cell_nodes[1](0, 3, 1) = 1.0;
00300     cell_nodes[1](0, 3, 2) = -1.0;
00301     cell_nodes[1](0, 4, 0) = -1.0;
00302     cell_nodes[1](0, 4, 1) = -1.0;
00303     cell_nodes[1](0, 4, 2) = 1.0;
00304     cell_nodes[1](0, 5, 0) = 1.0;
00305     cell_nodes[1](0, 5, 1) = -1.0;
00306     cell_nodes[1](0, 5, 2) = 1.0;
00307     cell_nodes[1](0, 6, 0) = 1.0;
00308     cell_nodes[1](0, 6, 1) = 1.0;
00309     cell_nodes[1](0, 6, 2) = 1.0;
00310     cell_nodes[1](0, 7, 0) = -1.0;
00311     cell_nodes[1](0, 7, 1) = 1.0;
00312     cell_nodes[1](0, 7, 2) = 1.0;
00313 
00314     std::stringstream mystream[2];
00315     mystream[0].str("\n>> Now testing basis on a generic parallelepiped ...\n");
00316     mystream[1].str("\n>> Now testing basis on the reference hex ...\n");
00317 
00318 
00319     for (int pcell = 0; pcell < 2; pcell++) {
00320       *outStream << mystream[pcell].str();
00321       FieldContainer<double> interp_points(1, numInterpPoints, cellDim);
00322       CellTools<double>::mapToPhysicalFrame(interp_points, interp_points_ref, cell_nodes[pcell], cell);
00323       interp_points.resize(numInterpPoints, cellDim);
00324 
00325       for (int x_order=0; x_order <= max_order; x_order++) {
00326         int max_y_order = max_order;
00327         if (pcell == 0) {
00328           max_y_order -= x_order;
00329         }
00330         for (int y_order=0; y_order <= max_y_order; y_order++) {
00331           int max_z_order = max_order;
00332           if (pcell == 0) {
00333             max_z_order -= x_order;
00334             max_z_order -= y_order;
00335           }
00336           for (int z_order=0; z_order <= max_z_order; z_order++) {
00337 
00338             // evaluate exact solution
00339             FieldContainer<double> exact_solution(1, numInterpPoints);
00340             u_exact(exact_solution, interp_points, x_order, y_order, z_order);
00341 
00342             int basis_order = max_order;
00343 
00344             // set test tolerance;
00345             double zero = basis_order*basis_order*basis_order*100*INTREPID_TOL;
00346 
00347             //create basis
00348             FieldContainer<double> pts(PointTools::getLatticeSize(line,basis_order),1);
00349             PointTools::getLattice<double,FieldContainer<double> >(pts,line,basis_order);
00350             
00351             Teuchos::RCP<Basis<double,FieldContainer<double> > > basis =
00352               Teuchos::rcp(new Basis_HGRAD_HEX_Cn_FEM<double,FieldContainer<double> >( basis_order, POINTTYPE_SPECTRAL ) );
00353             int numFields = basis->getCardinality();
00354 
00355             // create cubatures
00356             Teuchos::RCP<Cubature<double> > cellCub = cubFactory.create(cell, 2*basis_order);
00357             Teuchos::RCP<Cubature<double> > sideCub = cubFactory.create(side, 2*basis_order);
00358             int numCubPointsCell = cellCub->getNumPoints();
00359             int numCubPointsSide = sideCub->getNumPoints();
00360 
00361             /* Computational arrays. */
00362             /* Section 1: Related to parent cell integration. */
00363             FieldContainer<double> cub_points_cell(numCubPointsCell, cellDim);
00364             FieldContainer<double> cub_points_cell_physical(1, numCubPointsCell, cellDim);
00365             FieldContainer<double> cub_weights_cell(numCubPointsCell);
00366             FieldContainer<double> jacobian_cell(1, numCubPointsCell, cellDim, cellDim);
00367             FieldContainer<double> jacobian_inv_cell(1, numCubPointsCell, cellDim, cellDim);
00368             FieldContainer<double> jacobian_det_cell(1, numCubPointsCell);
00369             FieldContainer<double> weighted_measure_cell(1, numCubPointsCell);
00370 
00371             FieldContainer<double> value_of_basis_at_cub_points_cell(numFields, numCubPointsCell);
00372             FieldContainer<double> transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00373             FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell);
00374             FieldContainer<double> grad_of_basis_at_cub_points_cell(numFields, numCubPointsCell, cellDim);
00375             FieldContainer<double> transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00376             FieldContainer<double> weighted_transformed_grad_of_basis_at_cub_points_cell(1, numFields, numCubPointsCell, cellDim);
00377             FieldContainer<double> fe_matrix(1, numFields, numFields);
00378 
00379             FieldContainer<double> rhs_at_cub_points_cell_physical(1, numCubPointsCell);
00380             FieldContainer<double> rhs_and_soln_vector(1, numFields);
00381 
00382             /* Section 2: Related to subcell (side) integration. */
00383             FieldContainer<double> cub_points_side(numCubPointsSide, sideDim);
00384             FieldContainer<double> cub_weights_side(numCubPointsSide);
00385             FieldContainer<double> cub_points_side_refcell(numCubPointsSide, cellDim);
00386             FieldContainer<double> cub_points_side_physical(1, numCubPointsSide, cellDim);
00387             FieldContainer<double> jacobian_side_refcell(1, numCubPointsSide, cellDim, cellDim);
00388             FieldContainer<double> jacobian_det_side_refcell(1, numCubPointsSide);
00389             FieldContainer<double> weighted_measure_side_refcell(1, numCubPointsSide);
00390 
00391             FieldContainer<double> value_of_basis_at_cub_points_side_refcell(numFields, numCubPointsSide);
00392             FieldContainer<double> transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00393             FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points_side_refcell(1, numFields, numCubPointsSide);
00394             FieldContainer<double> neumann_data_at_cub_points_side_physical(1, numCubPointsSide);
00395             FieldContainer<double> neumann_fields_per_side(1, numFields);
00396 
00397             /* Section 3: Related to global interpolant. */
00398             FieldContainer<double> value_of_basis_at_interp_points_ref(numFields, numInterpPoints);
00399             FieldContainer<double> transformed_value_of_basis_at_interp_points_ref(1, numFields, numInterpPoints);
00400             FieldContainer<double> interpolant(1, numInterpPoints);
00401 
00402             FieldContainer<int> ipiv(numFields);
00403 
00404 
00405 
00406             /******************* START COMPUTATION ***********************/
00407 
00408             // get cubature points and weights
00409             cellCub->getCubature(cub_points_cell, cub_weights_cell);
00410 
00411             // compute geometric cell information
00412             CellTools<double>::setJacobian(jacobian_cell, cub_points_cell, cell_nodes[pcell], cell);
00413             CellTools<double>::setJacobianInv(jacobian_inv_cell, jacobian_cell);
00414             CellTools<double>::setJacobianDet(jacobian_det_cell, jacobian_cell);
00415 
00416             // compute weighted measure
00417             FunctionSpaceTools::computeCellMeasure<double>(weighted_measure_cell, jacobian_det_cell, cub_weights_cell);
00418 
00420             // Computing mass matrices:
00421             // tabulate values of basis functions at (reference) cubature points
00422             basis->getValues(value_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_VALUE);
00423 
00424             // transform values of basis functions 
00425             FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_cell,
00426                                                             value_of_basis_at_cub_points_cell);
00427 
00428             // multiply with weighted measure
00429             FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_cell,
00430                                                         weighted_measure_cell,
00431                                                         transformed_value_of_basis_at_cub_points_cell);
00432 
00433             // compute mass matrices
00434             FunctionSpaceTools::integrate<double>(fe_matrix,
00435                                                   transformed_value_of_basis_at_cub_points_cell,
00436                                                   weighted_transformed_value_of_basis_at_cub_points_cell,
00437                                                   COMP_BLAS);
00439 
00441             // Computing stiffness matrices:
00442             // tabulate gradients of basis functions at (reference) cubature points
00443             basis->getValues(grad_of_basis_at_cub_points_cell, cub_points_cell, OPERATOR_GRAD);
00444 
00445             // transform gradients of basis functions 
00446             FunctionSpaceTools::HGRADtransformGRAD<double>(transformed_grad_of_basis_at_cub_points_cell,
00447                                                            jacobian_inv_cell,
00448                                                            grad_of_basis_at_cub_points_cell);
00449 
00450             // multiply with weighted measure
00451             FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_grad_of_basis_at_cub_points_cell,
00452                                                         weighted_measure_cell,
00453                                                         transformed_grad_of_basis_at_cub_points_cell);
00454 
00455             // compute stiffness matrices and sum into fe_matrix
00456             FunctionSpaceTools::integrate<double>(fe_matrix,
00457                                                   transformed_grad_of_basis_at_cub_points_cell,
00458                                                   weighted_transformed_grad_of_basis_at_cub_points_cell,
00459                                                   COMP_BLAS,
00460                                                   true);
00462 
00464             // Computing RHS contributions:
00465             // map cell (reference) cubature points to physical space
00466             CellTools<double>::mapToPhysicalFrame(cub_points_cell_physical, cub_points_cell, cell_nodes[pcell], cell);
00467 
00468             // evaluate rhs function
00469             rhsFunc(rhs_at_cub_points_cell_physical, cub_points_cell_physical, x_order, y_order, z_order);
00470 
00471             // compute rhs
00472             FunctionSpaceTools::integrate<double>(rhs_and_soln_vector,
00473                                                   rhs_at_cub_points_cell_physical,
00474                                                   weighted_transformed_value_of_basis_at_cub_points_cell,
00475                                                   COMP_BLAS);
00476 
00477             // compute neumann b.c. contributions and adjust rhs
00478             sideCub->getCubature(cub_points_side, cub_weights_side);
00479             for (unsigned i=0; i<numSides; i++) {
00480               // compute geometric cell information
00481               CellTools<double>::mapToReferenceSubcell(cub_points_side_refcell, cub_points_side, sideDim, (int)i, cell);
00482               CellTools<double>::setJacobian(jacobian_side_refcell, cub_points_side_refcell, cell_nodes[pcell], cell);
00483               CellTools<double>::setJacobianDet(jacobian_det_side_refcell, jacobian_side_refcell);
00484 
00485               // compute weighted face measure
00486               FunctionSpaceTools::computeFaceMeasure<double>(weighted_measure_side_refcell,
00487                                                              jacobian_side_refcell,
00488                                                              cub_weights_side,
00489                                                              i,
00490                                                              cell);
00491 
00492               // tabulate values of basis functions at side cubature points, in the reference parent cell domain
00493               basis->getValues(value_of_basis_at_cub_points_side_refcell, cub_points_side_refcell, OPERATOR_VALUE);
00494               // transform 
00495               FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_cub_points_side_refcell,
00496                                                               value_of_basis_at_cub_points_side_refcell);
00497 
00498               // multiply with weighted measure
00499               FunctionSpaceTools::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00500                                                           weighted_measure_side_refcell,
00501                                                           transformed_value_of_basis_at_cub_points_side_refcell);
00502 
00503               // compute Neumann data
00504               // map side cubature points in reference parent cell domain to physical space
00505               CellTools<double>::mapToPhysicalFrame(cub_points_side_physical, cub_points_side_refcell, cell_nodes[pcell], cell);
00506               // now compute data
00507               neumann(neumann_data_at_cub_points_side_physical, cub_points_side_physical, jacobian_side_refcell,
00508                       cell, (int)i, x_order, y_order, z_order);
00509 
00510               FunctionSpaceTools::integrate<double>(neumann_fields_per_side,
00511                                                     neumann_data_at_cub_points_side_physical,
00512                                                     weighted_transformed_value_of_basis_at_cub_points_side_refcell,
00513                                                     COMP_BLAS);
00514 
00515               // adjust RHS
00516               RealSpaceTools<double>::add(rhs_and_soln_vector, neumann_fields_per_side);;
00517             }
00519 
00521             // Solution of linear system:
00522             int info = 0;
00523             Teuchos::LAPACK<int, double> solver;
00524             solver.GESV(numFields, 1, &fe_matrix[0], numFields, &ipiv(0), &rhs_and_soln_vector[0], numFields, &info);
00526 
00528             // Building interpolant:
00529             // evaluate basis at interpolation points
00530             basis->getValues(value_of_basis_at_interp_points_ref, interp_points_ref, OPERATOR_VALUE);
00531             // transform values of basis functions 
00532             FunctionSpaceTools::HGRADtransformVALUE<double>(transformed_value_of_basis_at_interp_points_ref,
00533                                                             value_of_basis_at_interp_points_ref);
00534             FunctionSpaceTools::evaluate<double>(interpolant, rhs_and_soln_vector, transformed_value_of_basis_at_interp_points_ref);
00536 
00537             /******************* END COMPUTATION ***********************/
00538         
00539             RealSpaceTools<double>::subtract(interpolant, exact_solution);
00540 
00541             *outStream << "\nRelative norm-2 error between exact solution polynomial of order ("
00542                        << x_order << ", " << y_order << ", " << z_order
00543                        << ") and finite element interpolant of order " << basis_order << ": "
00544                        << RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00545                           RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) << "\n";
00546 
00547             if (RealSpaceTools<double>::vectorNorm(&interpolant[0], interpolant.dimension(1), NORM_TWO) /
00548                 RealSpaceTools<double>::vectorNorm(&exact_solution[0], exact_solution.dimension(1), NORM_TWO) > zero) {
00549               *outStream << "\n\nPatch test failed for solution polynomial order ("
00550                          << x_order << ", " << y_order << ", " << z_order << ") and basis order " << basis_order << "\n\n";
00551               errorFlag++;
00552             }
00553           } // end for z_order
00554         } // end for y_order
00555       } // end for x_order
00556     } // end for pcell
00557 
00558   }
00559   // Catch unexpected errors
00560   catch (std::logic_error err) {
00561     *outStream << err.what() << "\n\n";
00562     errorFlag = -1000;
00563   };
00564 
00565   if (errorFlag != 0)
00566     std::cout << "End Result: TEST FAILED\n";
00567   else
00568     std::cout << "End Result: TEST PASSED\n";
00569 
00570   // reset format state of std::cout
00571   std::cout.copyfmt(oldFormatState);
00572 
00573   return errorFlag;
00574 }