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