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