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