SundanceProblemTesting.cpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                             Sundance
00005 //                 Copyright 2011 Sandia Corporation
00006 // 
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Kevin Long (kevin.long@ttu.edu)
00038 // 
00039 
00040 /* @HEADER@ */
00041 
00042 #include "SundanceProblemTesting.hpp"
00043 
00044 #include "SundanceMaximalCellFilter.hpp"
00045 #include "SundanceDimensionalCellFilter.hpp"
00046 #include "SundancePositionalCellPredicate.hpp"
00047 #include "SundanceMesh.hpp"
00048 #include "SundanceMeshSource.hpp"
00049 #include "SundanceMeshType.hpp"
00050 #include "SundanceBasicSimplicialMeshType.hpp"
00051 #include "SundancePartitionedLineMesher.hpp"
00052 #include "SundancePartitionedRectangleMesher.hpp"
00053 #include "PlayaEpetraVectorType.hpp"
00054 #include "PlayaLinearSolverBuilder.hpp"
00055 #include "SundanceGaussianQuadrature.hpp"
00056 #include "SundanceCoordExpr.hpp"
00057 #include "SundanceCellDiameterExpr.hpp"
00058 #include "SundanceLinearProblem.hpp"
00059 #include "SundanceIntegral.hpp"
00060 
00061 namespace Sundance
00062 {
00063 
00064   using std::setprecision;
00065   using std::setw;
00066 
00067 bool checkErrorNorms(
00068   const Mesh& mesh,
00069   const CellFilter& filter,
00070   const Expr& numSoln,
00071   const Expr& exactSoln,
00072   const QuadratureFamily& quad,
00073   double L2Tol,
00074   double H1SemiTol,
00075   double H1Tol)
00076 {
00077   Tabs tab;
00078   Tabs tab1;
00079   Expr df = numSoln - exactSoln;
00080 
00081   double L2Err = L2Norm(mesh, filter, df, quad);
00082   Out::root() << tab << "L2 norm check:" << std::endl
00083               << tab1 << setw(16) << setprecision(6) << L2Err 
00084               << " tol=" << setw(12) << L2Tol ;
00085   if (fabs(L2Err) > L2Tol) Out::root() << "   <==== FAIL!" << std::endl;
00086   else Out::root() << std::endl;
00087 
00088   double H1SemiErr = H1Seminorm(mesh, filter, df, quad);
00089   Out::root() << tab << "H1 seminorm check:" << std::endl
00090               << tab1 << setw(16) << setprecision(6) << H1SemiErr 
00091               << " tol=" << setw(12) << H1SemiTol ;
00092   if (fabs(H1SemiErr) > H1SemiTol) Out::root() << "   <==== FAIL!" << std::endl;
00093   else Out::root() << std::endl;
00094 
00095 
00096   double H1Err = H1Norm(mesh, filter, df, quad);
00097   Out::root() << tab << "H1 norm check:" << std::endl
00098               << tab1 << setw(16) << setprecision(6) << H1Err 
00099               << " tol=" << setw(12) << H1Tol;
00100   if (fabs(H1Err) > H1Tol) Out::root() << "   <==== FAIL!" << std::endl;
00101   else Out::root() << std::endl;
00102 
00103   return (fabs(L2Err) <= L2Tol) 
00104     && (fabs(H1SemiErr) <= H1SemiTol) 
00105     && (fabs(H1Err) <= H1Tol) ;
00106 }
00107 
00108 double fitPower(const Array<double>& h, const Array<double>& err)
00109 {
00110   Array<double> x(h.size());
00111   Array<double> y(h.size());
00112   double xBar = 0.0;
00113   double yBar = 0.0;
00114   for (int i=0; i<h.size(); i++)
00115   {
00116     x[i] = log(h[i]);
00117     y[i] = log(err[i]);
00118     xBar += x[i];
00119     yBar += y[i];
00120   }
00121 
00122   xBar /= h.size();
00123   yBar /= h.size();
00124   
00125   double u = 0.0;
00126   double v = 0.0;
00127   for (int i=0; i<h.size(); i++)
00128   {
00129     u += (x[i]-xBar)*(y[i]-yBar);
00130     v += pow(x[i]-xBar,2.0);
00131   }
00132 
00133   return u/v;
00134 }
00135 
00136 /* -------- LineDomain ----------------- */
00137 
00138 LineDomain::LineDomain(const Array<int>& nx)
00139   : a_(0.0), b_(1.0), nx_(nx), interior_(new MaximalCellFilter()),
00140     left_(), right_(), mesh_()
00141 {init();}
00142 
00143 
00144 LineDomain::LineDomain(double a, double b, const Array<int>& nx)
00145   : a_(a), b_(b), nx_(nx), interior_(new MaximalCellFilter()),
00146     left_(), right_(), mesh_()
00147 {init();}
00148 
00149 void LineDomain::init()
00150 {
00151   int np = MPIComm::world().getNProc();
00152   MeshType meshType = new BasicSimplicialMeshType();
00153   mesh_.resize(nx_.size());
00154   for (int i=0; i<nx_.size(); i++)
00155   {
00156     MeshSource mesher = new PartitionedLineMesher(a_, b_, np*nx_[i], meshType);
00157     mesh_[i] = mesher.getMesh();
00158   }
00159   
00160   CellFilter points = new DimensionalCellFilter(0);
00161   left_ = points.subset(new CoordinateValueCellPredicate(0,a_));
00162   right_ = points.subset(new CoordinateValueCellPredicate(0,b_));
00163 }
00164 
00165 /* -------- RectangleDomain ----------------- */
00166 
00167 RectangleDomain::RectangleDomain(const Array<int>& n)
00168   : ax_(0.0), bx_(1.0), nx_(n), 
00169     ay_(0.0), by_(1.0), ny_(n), 
00170     interior_(new MaximalCellFilter()),
00171     north_(), south_(), east_(), west_(),
00172     mesh_()
00173 {init();}
00174 
00175 RectangleDomain::RectangleDomain(
00176   double ax, double bx, const Array<int>& nx,
00177   double ay, double by, const Array<int>& ny
00178   )
00179   : ax_(ax), bx_(bx), nx_(nx), 
00180     ay_(ay), by_(by), ny_(ny), 
00181     interior_(new MaximalCellFilter()),
00182     north_(), south_(), east_(), west_(),
00183     mesh_()
00184 {init();}
00185 
00186 
00187 void RectangleDomain::init()
00188 {
00189   TEUCHOS_TEST_FOR_EXCEPT(nx_.size() != ny_.size());
00190   int np = MPIComm::world().getNProc();
00191   int npx = -1;
00192   int npy = -1;
00193 
00194   PartitionedRectangleMesher::balanceXY(np, &npx, &npy);
00195   TEUCHOS_TEST_FOR_EXCEPT(npx < 1);
00196   TEUCHOS_TEST_FOR_EXCEPT(npy < 1);
00197   TEUCHOS_TEST_FOR_EXCEPT(npx * npy != np);
00198 
00199   MeshType meshType = new BasicSimplicialMeshType();
00200   mesh_.resize(nx_.size());
00201 
00202   for (int i=0; i<nx_.size(); i++)
00203   {
00204     MeshSource mesher 
00205       = new PartitionedRectangleMesher(
00206         ax_, bx_, nx_[i], npx,
00207         ay_, by_, ny_[i], npy,
00208         meshType);
00209     mesh_[i] = mesher.getMesh();
00210   }
00211   
00212   CellFilter edges = new DimensionalCellFilter(1);
00213   north_ = edges.subset(new CoordinateValueCellPredicate(1,by_));
00214   south_ = edges.subset(new CoordinateValueCellPredicate(1,ay_));
00215   west_ = edges.subset(new CoordinateValueCellPredicate(0,ax_));
00216   east_ = edges.subset(new CoordinateValueCellPredicate(0,bx_));
00217 }
00218 
00219 
00220 
00221 
00222 
00223 /* -------- norm calculators ----------------- */
00224 
00225 Array<double> L2NormCalculator::computeNorms(
00226   const ForwardProblemTestBase* prob,
00227   int meshIndex,
00228   const Expr& numSoln, const Expr& exactSoln) const
00229 {
00230   Expr errFunc = (numSoln - exactSoln).flatten();
00231 
00232   Mesh mesh = prob->getMesh(meshIndex);
00233   CellFilter interior = prob->interior();
00234 
00235   Array<int> p = prob->pExpected();
00236   TEUCHOS_TEST_FOR_EXCEPTION(p.size() != errFunc.size(),
00237     std::runtime_error,
00238     "size mismatch between array of expected orders (p=" << p << ") and "
00239     "array of solutions: " << errFunc.size());
00240  
00241   Array<double> rtn(p.size());
00242   for (int i=0; i<errFunc.size(); i++)
00243   {
00244     QuadratureFamily quad = new GaussianQuadrature(2*p[i]);
00245     double L2Err = L2Norm(mesh, interior, errFunc[i], quad);
00246     rtn[i] = L2Err;
00247   }  
00248 
00249   return rtn;
00250 }
00251   
00252 
00253 Array<LPTestSpec> LPTestBase::specs() const
00254 {
00255   return tuple(
00256     LPTestSpec("amesos.xml", 1.0e-10, makeSet<int>(1)),
00257     LPTestSpec("aztec-ifpack.xml", 1.0e-10),
00258     LPTestSpec("aztec-ml.xml", 1.0e-10),
00259     LPTestSpec("belos-ifpack.xml", 1.0e-8),
00260     LPTestSpec("belos-ml.xml", 1.0e-10)
00261     );
00262 }
00263 
00264 /** \relates LPTestSpec */
00265 std::ostream& operator<<(std::ostream& os, const LPTestSpec& spec)
00266 {
00267   os << "LPTestSpec(tol=" << spec.tol() << ", solver=" << spec.solverFile()
00268      << std::endl;
00269   return os;
00270 }
00271 
00272 
00273 
00274 
00275 bool LPTestSuite::run() const 
00276 {
00277   int np = MPIComm::world().getNProc();
00278 
00279   bool allOK = true;
00280 
00281   for (int i=0; i<tests_.size(); i++)
00282   {
00283     Tabs tab(0);
00284     Array<LPTestSpec> specs = tests_[i]->specs();
00285     for (int j=0; j<specs.size(); j++)
00286     {
00287       if (specs[j].nProcIsAllowed(np))
00288       {
00289         Out::root() << std::endl;
00290         Out::root() << std::endl;
00291         Out::root() << std::endl;
00292         Out::root() << tab 
00293                     << "-------------------------------------" 
00294           "-------------------------------------" 
00295                     << std::endl;
00296         Out::root() << tab << "running test " << tests_[i]->name()
00297                     << " with spec " << specs[j] << std::endl;
00298         Out::root() << tab 
00299                     << "-------------------------------------" 
00300           "-------------------------------------" 
00301                     << std::endl;
00302 
00303         std::string solverFile = specs[j].solverFile();
00304         double tol = specs[j].tol();
00305         bool pass = tests_[i]->run(solverFile, tol);
00306         allOK = pass && allOK;
00307       }
00308       else
00309       {
00310         Out::root() << tab << "skipping test " << tests_[i]->name()
00311                     << " with spec=" << specs[j] << std::endl;
00312       }
00313     }
00314     Out::root() << std::endl;
00315     Out::root() << std::endl;
00316   }
00317   return allOK;
00318 }
00319 
00320 
00321 LPTestSuite::LPTestSuite()
00322   : tests_(), testSpecs_() {}
00323 
00324 
00325 void LPTestSuite::registerTest(const RCP<LPTestBase>& test) 
00326 {
00327   tests_.append(test);
00328 }
00329 
00330 VectorType<double> ForwardProblemTestBase::vecType() const
00331 {
00332   return new EpetraVectorType();
00333 }
00334 
00335 Expr ForwardProblemTestBase::coord(int d) const 
00336 {
00337   TEUCHOS_TEST_FOR_EXCEPT(d<0 || d>2);
00338   return new CoordExpr(d);
00339 }
00340 
00341 double ForwardProblemTestBase::cellSize(int i) const 
00342 {
00343   Expr h = new CellDiameterExpr();
00344   Expr hExpr = Integral(interior(), h, new GaussianQuadrature(1));
00345   Expr AExpr = Integral(interior(), 1.0, new GaussianQuadrature(1));
00346   
00347   double area = evaluateIntegral(getMesh(i), AExpr);
00348   double hMean = evaluateIntegral(getMesh(i), hExpr)/area;
00349 
00350   return hMean;
00351 }
00352 
00353 RCP<ErrNormCalculatorBase> ForwardProblemTestBase::normCalculator() const
00354 {
00355   return rcp(new L2NormCalculator());
00356 
00357 }
00358 bool ForwardProblemTestBase::run(const std::string& solverFile, 
00359   double tol) const
00360 {
00361   if (numMeshes()==1) return runSingleTest(solverFile, tol);
00362   else return runTestSequence(solverFile, tol);
00363 }
00364 
00365 
00366 bool ForwardProblemTestBase::runTestSequence(const std::string& solverFile, 
00367   const double& pTol) const
00368 {
00369   Tabs tab(0);
00370   Out::root() << tab << "Running test sequence " << name() << std::endl;
00371 
00372   LinearSolver<double> solver 
00373     = LinearSolverBuilder::createSolver(solverFile);
00374 
00375   Expr exact = exactSoln();
00376 
00377   Array<double> hList(numMeshes());
00378   Array<Array<double> > errList(numMeshes());
00379 
00380   for (int i=0; i<numMeshes(); i++)
00381   {
00382     Tabs tab1;
00383     Out::root() << tab1 << "running on mesh #" << i << std::endl;
00384     Expr soln;
00385     bool solveOK = solve(getMesh(i), solver, soln);
00386     if (!solveOK) return false;
00387 
00388     postRunCallback(i, getMesh(i), solverFile, soln);
00389 
00390     double h = cellSize(i);
00391     Array<double> err = normCalculator()->computeNorms(
00392       this, i, soln, exact);
00393     hList[i] = h;
00394     errList[i] = err;
00395   }
00396 
00397   bool allPass = true;
00398   Out::root() << tab << "Error results: " << std::endl;
00399   for (int f=0; f < pExpected().size(); f++)
00400   {
00401     Tabs tab1;
00402     Out::root() << tab1 << "Variable #" << f << std::endl; 
00403     Out::root() << std::endl << tab1 << "Error norm versus h" << std::endl;
00404     Array<double> err;
00405     for (int i=0; i<numMeshes(); i++)
00406     {
00407       Tabs tab2;
00408       Out::root() << tab2 << setw(20) << setprecision(5) << hList[i] 
00409                   << " " << setw(20) << setprecision(5) << errList[i][f] 
00410                   << std::endl;
00411       err.append(errList[i][f]); 
00412     }
00413     
00414     double p = fitPower(hList, err);
00415     Out::root() << tab << "Measured exponent: " << p << std::endl;
00416     Out::root() << tab << "Expected exponent: " << pExpected()[f] << std::endl;
00417     double pErr = fabs(p - pExpected()[f]) ;
00418     Out::root() << tab << "Difference: " << setw(12) << pErr
00419                 << " Tolerance: " << setw(12) << pTol;
00420     bool pass = pErr <= pTol;
00421     if (!pass) Out::root() << "  <==== FAIL!";
00422     Out::root() << std::endl;
00423     allPass = allPass && pass;
00424   }
00425 
00426   return allPass;
00427 }
00428 
00429 
00430 
00431 bool ForwardProblemTestBase::runSingleTest(const std::string& solverFile, 
00432   const double& tol) const
00433 {
00434   Tabs tab(0);
00435 
00436   LinearSolver<double> solver 
00437     = LinearSolverBuilder::createSolver(solverFile);
00438 
00439   Expr exact = exactSoln();
00440 
00441   Expr soln;
00442   bool solveOK = solve(getMesh(0), solver, soln);
00443   if (!solveOK) return false;
00444   
00445   Array<double> err = normCalculator()->computeNorms(this, 0, 
00446     soln, exact);
00447   
00448   for (int i=0; i<err.size(); i++)
00449   {
00450     Out::root() << tab << setw(20) << setprecision(5) << err[i] 
00451                 << " " << setw(20) << setprecision(5) << tol 
00452               << std::endl;
00453   }
00454   return err[0] <= tol;
00455 }
00456 
00457 bool LPTestBase::solve(const Mesh& mesh,
00458   const LinearSolver<double>& solver,
00459   Expr& soln) const
00460 {
00461   Tabs tab(0);
00462   LinearProblem::solveFailureIsFatal() = false;
00463   SolverState<double> state = prob(mesh).solve(solver, soln);
00464 
00465   bool converged = (state.finalState() == SolveConverged) ;
00466   if (!converged)
00467   {
00468     Out::root() << tab << "solve failed to converge!" << std::endl;
00469     Tabs tab1;
00470     Out::root() << tab1 << "state = " << state << std::endl;
00471   }
00472   return converged;
00473 }
00474 
00475 LP1DTestBase::LP1DTestBase(const Array<int>& nx)
00476   : domain_(nx) {}
00477 
00478 LP1DTestBase::LP1DTestBase(double a, double b, const Array<int>& nx)
00479   : domain_(a, b, nx) {}
00480 
00481 LPRectTestBase::LPRectTestBase(const Array<int>& n)
00482   : domain_(n) {}
00483 }

Site Contact