SundanceProblemTesting.cpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                              Sundance
00005 //                 Copyright (2005) Sandia Corporation
00006 // 
00007 // Copyright (year first published) Sandia Corporation.  Under the terms 
00008 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government 
00009 // retains certain rights in this software.
00010 // 
00011 // This library is free software; you can redistribute it and/or modify
00012 // it under the terms of the GNU Lesser General Public License as
00013 // published by the Free Software Foundation; either version 2.1 of the
00014 // License, or (at your option) any later version.
00015 //  
00016 // This library is distributed in the hope that it will be useful, but
00017 // WITHOUT ANY WARRANTY; without even the implied warranty of
00018 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00019 // Lesser General Public License for more details.
00020 //                                                                                 
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License along with this library; if not, write to the Free Software
00023 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00024 // USA                                                                                
00025 // Questions? Contact Kevin Long (krlong@sandia.gov), 
00026 // Sandia National Laboratories, Livermore, California, USA
00027 // 
00028 // ************************************************************************
00029 /* @HEADER@ */
00030 
00031 #include "SundanceProblemTesting.hpp"
00032 
00033 #include "SundanceMaximalCellFilter.hpp"
00034 #include "SundanceDimensionalCellFilter.hpp"
00035 #include "SundancePositionalCellPredicate.hpp"
00036 #include "SundanceMesh.hpp"
00037 #include "SundanceMeshSource.hpp"
00038 #include "SundanceMeshType.hpp"
00039 #include "SundanceBasicSimplicialMeshType.hpp"
00040 #include "SundancePartitionedLineMesher.hpp"
00041 #include "SundancePartitionedRectangleMesher.hpp"
00042 #include "TSFEpetraVectorType.hpp"
00043 #include "TSFLinearSolverBuilder.hpp"
00044 #include "SundanceGaussianQuadrature.hpp"
00045 #include "SundanceCoordExpr.hpp"
00046 #include "SundanceCellDiameterExpr.hpp"
00047 #include "SundanceLinearProblem.hpp"
00048 #include "SundanceIntegral.hpp"
00049 
00050 namespace Sundance
00051 {
00052 
00053 bool checkErrorNorms(
00054   const Mesh& mesh,
00055   const CellFilter& filter,
00056   const Expr& numSoln,
00057   const Expr& exactSoln,
00058   const QuadratureFamily& quad,
00059   double L2Tol,
00060   double H1SemiTol,
00061   double H1Tol)
00062 {
00063   Tabs tab;
00064   Tabs tab1;
00065   Expr df = numSoln - exactSoln;
00066 
00067   double L2Err = L2Norm(mesh, filter, df, quad);
00068   Out::root() << tab << "L2 norm check:" << std::endl
00069               << tab1 << setw(16) << setprecision(6) << L2Err 
00070               << " tol=" << setw(12) << L2Tol ;
00071   if (fabs(L2Err) > L2Tol) Out::root() << "   <==== FAIL!" << std::endl;
00072   else Out::root() << std::endl;
00073 
00074   double H1SemiErr = H1Seminorm(mesh, filter, df, quad);
00075   Out::root() << tab << "H1 seminorm check:" << std::endl
00076               << tab1 << setw(16) << setprecision(6) << H1SemiErr 
00077               << " tol=" << setw(12) << H1SemiTol ;
00078   if (fabs(H1SemiErr) > H1SemiTol) Out::root() << "   <==== FAIL!" << std::endl;
00079   else Out::root() << std::endl;
00080 
00081 
00082   double H1Err = H1Norm(mesh, filter, df, quad);
00083   Out::root() << tab << "H1 norm check:" << std::endl
00084               << tab1 << setw(16) << setprecision(6) << H1Err 
00085               << " tol=" << setw(12) << H1Tol;
00086   if (fabs(H1Err) > H1Tol) Out::root() << "   <==== FAIL!" << std::endl;
00087   else Out::root() << std::endl;
00088 
00089   return (fabs(L2Err) <= L2Tol) 
00090     && (fabs(H1SemiErr) <= H1SemiTol) 
00091     && (fabs(H1Err) <= H1Tol) ;
00092 }
00093 
00094 double fitPower(const Array<double>& h, const Array<double>& err)
00095 {
00096   Array<double> x(h.size());
00097   Array<double> y(h.size());
00098   double xBar = 0.0;
00099   double yBar = 0.0;
00100   for (int i=0; i<h.size(); i++)
00101   {
00102     x[i] = log(h[i]);
00103     y[i] = log(err[i]);
00104     xBar += x[i];
00105     yBar += y[i];
00106   }
00107 
00108   xBar /= h.size();
00109   yBar /= h.size();
00110   
00111   double u = 0.0;
00112   double v = 0.0;
00113   for (int i=0; i<h.size(); i++)
00114   {
00115     u += (x[i]-xBar)*(y[i]-yBar);
00116     v += pow(x[i]-xBar,2.0);
00117   }
00118 
00119   return u/v;
00120 }
00121 
00122 /* -------- LineDomain ----------------- */
00123 
00124 LineDomain::LineDomain(const Array<int>& nx)
00125   : a_(0.0), b_(1.0), nx_(nx), interior_(new MaximalCellFilter()),
00126     left_(), right_(), mesh_()
00127 {init();}
00128 
00129 
00130 LineDomain::LineDomain(double a, double b, const Array<int>& nx)
00131   : a_(a), b_(b), nx_(nx), interior_(new MaximalCellFilter()),
00132     left_(), right_(), mesh_()
00133 {init();}
00134 
00135 void LineDomain::init()
00136 {
00137   int np = MPIComm::world().getNProc();
00138   MeshType meshType = new BasicSimplicialMeshType();
00139   mesh_.resize(nx_.size());
00140   for (int i=0; i<nx_.size(); i++)
00141   {
00142     MeshSource mesher = new PartitionedLineMesher(a_, b_, np*nx_[i], meshType);
00143     mesh_[i] = mesher.getMesh();
00144   }
00145   
00146   CellFilter points = new DimensionalCellFilter(0);
00147   left_ = points.subset(new CoordinateValueCellPredicate(0,a_));
00148   right_ = points.subset(new CoordinateValueCellPredicate(0,b_));
00149 }
00150 
00151 
00152 
00153 Array<double> L2NormCalculator::computeNorms(
00154   const ForwardProblemTestBase* prob,
00155   int meshIndex,
00156   const Expr& numSoln, const Expr& exactSoln) const
00157 {
00158   Expr errFunc = (numSoln - exactSoln).flatten();
00159 
00160   Mesh mesh = prob->getMesh(meshIndex);
00161   CellFilter interior = prob->interior();
00162 
00163   Array<int> p = prob->pExpected();
00164   TEST_FOR_EXCEPTION(p.size() != errFunc.size(),
00165     RuntimeError,
00166     "size mismatch between array of expected orders (p=" << p << ") and "
00167     "array of solutions: " << errFunc.size());
00168  
00169   Array<double> rtn(p.size());
00170   for (int i=0; i<errFunc.size(); i++)
00171   {
00172     QuadratureFamily quad = new GaussianQuadrature(2*p[i]);
00173     double L2Err = L2Norm(mesh, interior, errFunc[i], quad);
00174     rtn[i] = L2Err;
00175   }  
00176 
00177   return rtn;
00178 }
00179   
00180 
00181 Array<LPTestSpec> LPTestBase::specs() const
00182 {
00183   return tuple(
00184     LPTestSpec("amesos.xml", 1.0e-10, makeSet<int>(1)),
00185     LPTestSpec("aztec-ifpack.xml", 1.0e-10),
00186     LPTestSpec("aztec-ml.xml", 1.0e-10),
00187     LPTestSpec("belos-ifpack.xml", 1.0e-8),
00188     LPTestSpec("belos-ml.xml", 1.0e-10)
00189     );
00190 }
00191 
00192 /** \relates LPTestSpec */
00193 std::ostream& operator<<(std::ostream& os, const LPTestSpec& spec)
00194 {
00195   os << "LPTestSpec(tol=" << spec.tol() << ", solver=" << spec.solverFile()
00196      << std::endl;
00197   return os;
00198 }
00199 
00200 
00201 
00202 
00203 bool LPTestSuite::run() const 
00204 {
00205   int np = MPIComm::world().getNProc();
00206 
00207   bool allOK = true;
00208 
00209   for (int i=0; i<tests_.size(); i++)
00210   {
00211     Tabs tab(0);
00212     Array<LPTestSpec> specs = tests_[i]->specs();
00213     for (int j=0; j<specs.size(); j++)
00214     {
00215       if (specs[j].nProcIsAllowed(np))
00216       {
00217         Out::root() << std::endl;
00218         Out::root() << std::endl;
00219         Out::root() << std::endl;
00220         Out::root() << tab 
00221                     << "-------------------------------------" 
00222           "-------------------------------------" 
00223                     << std::endl;
00224         Out::root() << tab << "running test " << tests_[i]->name()
00225                     << " with spec " << specs[j] << std::endl;
00226         Out::root() << tab 
00227                     << "-------------------------------------" 
00228           "-------------------------------------" 
00229                     << std::endl;
00230 
00231         std::string solverFile = specs[j].solverFile();
00232         double tol = specs[j].tol();
00233         bool pass = tests_[i]->run(solverFile, tol);
00234         allOK = pass && allOK;
00235       }
00236       else
00237       {
00238         Out::root() << tab << "skipping test " << tests_[i]->name()
00239                     << " with spec=" << specs[j] << std::endl;
00240       }
00241     }
00242     Out::root() << std::endl;
00243     Out::root() << std::endl;
00244   }
00245   return allOK;
00246 }
00247 
00248 
00249 LPTestSuite::LPTestSuite()
00250   : tests_(), testSpecs_() {}
00251 
00252 
00253 void LPTestSuite::registerTest(const RCP<LPTestBase>& test) 
00254 {
00255   tests_.append(test);
00256 }
00257 
00258 VectorType<double> ForwardProblemTestBase::vecType() const
00259 {
00260   return new EpetraVectorType();
00261 }
00262 
00263 Expr ForwardProblemTestBase::coord(int d) const 
00264 {
00265   TEST_FOR_EXCEPT(d<0 || d>2);
00266   return new CoordExpr(d);
00267 }
00268 
00269 double ForwardProblemTestBase::cellSize(int i) const 
00270 {
00271   Expr h = new CellDiameterExpr();
00272   Expr hExpr = Integral(interior(), h, new GaussianQuadrature(1));
00273   Expr AExpr = Integral(interior(), 1.0, new GaussianQuadrature(1));
00274   
00275   double area = evaluateIntegral(getMesh(i), AExpr);
00276   double hMean = evaluateIntegral(getMesh(i), hExpr)/area;
00277 
00278   return hMean;
00279 }
00280 
00281 RCP<ErrNormCalculatorBase> ForwardProblemTestBase::normCalculator() const
00282 {
00283   return rcp(new L2NormCalculator());
00284 
00285 }
00286 bool ForwardProblemTestBase::run(const std::string& solverFile, 
00287   double tol) const
00288 {
00289   if (numMeshes()==1) return runSingleTest(solverFile, tol);
00290   else return runTestSequence(solverFile, tol);
00291 }
00292 
00293 
00294 bool ForwardProblemTestBase::runTestSequence(const std::string& solverFile, 
00295   const double& pTol) const
00296 {
00297   Tabs tab(0);
00298   Out::root() << tab << "Running test sequence " << name() << std::endl;
00299 
00300   LinearSolver<double> solver 
00301     = LinearSolverBuilder::createSolver(solverFile);
00302 
00303   Expr exact = exactSoln();
00304 
00305   Array<double> hList(numMeshes());
00306   Array<Array<double> > errList(numMeshes());
00307 
00308   for (int i=0; i<numMeshes(); i++)
00309   {
00310     Tabs tab1;
00311     Out::root() << tab1 << "running on mesh #" << i << std::endl;
00312     Expr soln;
00313     bool solveOK = solve(getMesh(i), solver, soln);
00314     if (!solveOK) return false;
00315 
00316     double h = cellSize(i);
00317     Array<double> err = normCalculator()->computeNorms(
00318       this, i, soln, exact);
00319     hList[i] = h;
00320     errList[i] = err;
00321   }
00322 
00323   bool allPass = true;
00324   Out::root() << tab << "Error results: " << std::endl;
00325   for (int f=0; f < pExpected().size(); f++)
00326   {
00327     Tabs tab1;
00328     Out::root() << tab1 << "Variable #" << f << std::endl; 
00329     Out::root() << std::endl << tab1 << "Error norm versus h" << std::endl;
00330     Array<double> err;
00331     for (int i=0; i<numMeshes(); i++)
00332     {
00333       Tabs tab2;
00334       Out::root() << tab2 << setw(20) << setprecision(5) << hList[i] 
00335                   << " " << setw(20) << setprecision(5) << errList[i][f] 
00336                   << std::endl;
00337       err.append(errList[i][f]); 
00338     }
00339     
00340     double p = fitPower(hList, err);
00341     Out::root() << tab << "Measured exponent: " << p << std::endl;
00342     Out::root() << tab << "Expected exponent: " << pExpected()[f] << std::endl;
00343     double pErr = fabs(p - pExpected()[f]) ;
00344     Out::root() << tab << "Difference: " << setw(12) << pErr
00345                 << " Tolerance: " << setw(12) << pTol;
00346     bool pass = pErr <= pTol;
00347     if (!pass) Out::root() << "  <==== FAIL!";
00348     Out::root() << std::endl;
00349     allPass = allPass && pass;
00350   }
00351 
00352   return allPass;
00353 }
00354 
00355 
00356 
00357 bool ForwardProblemTestBase::runSingleTest(const std::string& solverFile, 
00358   const double& tol) const
00359 {
00360   Tabs tab(0);
00361 
00362   LinearSolver<double> solver 
00363     = LinearSolverBuilder::createSolver(solverFile);
00364 
00365   Expr exact = exactSoln();
00366 
00367   Expr soln;
00368   bool solveOK = solve(getMesh(0), solver, soln);
00369   if (!solveOK) return false;
00370   
00371   Array<double> err = normCalculator()->computeNorms(this, 0, 
00372     soln, exact);
00373   
00374   for (int i=0; i<err.size(); i++)
00375   {
00376     Out::root() << tab << setw(20) << setprecision(5) << err[i] 
00377                 << " " << setw(20) << setprecision(5) << tol 
00378               << std::endl;
00379   }
00380   return err[0] <= tol;
00381 }
00382 
00383 bool LPTestBase::solve(const Mesh& mesh,
00384   const LinearSolver<double>& solver,
00385   Expr& soln) const
00386 {
00387   Tabs tab(0);
00388   LinearProblem::solveFailureIsFatal() = false;
00389   SolverState<double> state = prob(mesh).solve(solver, soln);
00390 
00391   bool converged = (state.finalState() == SolveConverged) ;
00392   if (!converged)
00393   {
00394     Out::root() << tab << "solve failed to converge!" << std::endl;
00395     Tabs tab1;
00396     Out::root() << tab1 << "state = " << state << std::endl;
00397   }
00398   return converged;
00399 }
00400 
00401 LP1DTestBase::LP1DTestBase(const Array<int>& nx)
00402   : domain_(nx) {}
00403 
00404 LP1DTestBase::LP1DTestBase(double a, double b, const Array<int>& nx)
00405   : domain_(a, b, nx) {}
00406 
00407 }

Site Contact