LPTests1D.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 "Sundance.hpp"
00043 #include "SundanceCubicHermite.hpp"
00044 #include "SundanceProblemTesting.hpp"
00045 
00046 
00047 /** 
00048  * This class sets up a test that checks convergence rate for solution
00049  * of the Poisson equation in 1D with Lagrange(2) basis functions.
00050  * The test problem is 
00051  * \f[ u'' = -\frac{1}{4}\pi^2 \sin(\pi x/2) \f]
00052  * with boundary conditions
00053  * \f[ u(0)=1, \;\;\; u'(1)=0. \f]
00054  * The solution is \f$u(x)=\sin(\pi x/2)\f$.
00055  *
00056  * 
00057  */
00058 class SimplePoisson1DTest : public LP1DTestBase
00059 {
00060 public:
00061   /** 
00062    * Construct the test object. The constructor calls the base class
00063    * constructor with a list of mesh sizes to be used. 
00064    */
00065   SimplePoisson1DTest(const Array<int>& nx) 
00066     : LP1DTestBase(nx) {}
00067 
00068   /** Return a descriptive name for the test */
00069   std::string name() const {return "SimplePoisson1D";}
00070 
00071   /** Return the expected order of accuracy for the solutions. Because a
00072    * problem might have multiple field variables, the expected orders of 
00073    * accuracy for the different variables are returned as an array.  */
00074   Array<int> pExpected() const {return tuple<int>(3);}
00075 
00076   /** 
00077    * Return an expression for the exact solution to the test problem.
00078    */
00079   Expr exactSoln() const 
00080     {
00081       Expr x = coord(0);
00082       const double pi = 4.0*atan(1.0);
00083 
00084       return sin(pi*x/2.0);
00085     }
00086 
00087   /** 
00088    * Return a LP for the specified mesh. This function implements the
00089    * prob() pure virtual function of class LPTestBase.
00090    */
00091   LinearProblem prob(const Mesh& mesh) const
00092     {
00093       const double pi = 4.0*atan(1.0);
00094       CellFilter left = domain().left();
00095 
00096       Expr u = new UnknownFunction(new Lagrange(2), "u");
00097       Expr v = new TestFunction(new Lagrange(2), "v");
00098       Expr dx = gradient(1);
00099       Expr x = coord(0);
00100 
00101       QuadratureFamily quad = new GaussianQuadrature(4);
00102 
00103       Expr eqn = Integral(interior(), (dx*v)*(dx*u), quad)
00104         + Integral(interior(), -0.25*pi*pi*sin(pi*x/2.0)*v, quad);
00105 
00106       Expr bc = EssentialBC(left, v*u, quad);
00107       
00108       return LinearProblem(mesh, eqn, bc, v, u, vecType());
00109     }
00110 
00111   /**
00112    * Specify which linear solvers are to be used.
00113    */
00114   Array<LPTestSpec> specs() const
00115     {
00116       /* Use the standard solvers minus belos-ifpack, which has convergence
00117        * troubles on the finest mesh */
00118       double tol = 0.05;
00119       return tuple(
00120         LPTestSpec("amesos.xml", tol, makeSet<int>(1)),
00121         LPTestSpec("aztec-ifpack.xml", tol),
00122         LPTestSpec("aztec-ml.xml", tol),
00123         LPTestSpec("belos-ml.xml", tol),
00124         LPTestSpec("bicgstab.xml", tol)
00125         );
00126     }
00127 };
00128 
00129 
00130 
00131 
00132 /** Solves the Poisson equation in 1D with user-specified basis functions */
00133 class Poisson1DTest : public LP1DTestBase
00134 {
00135 public:
00136   /** */
00137   Poisson1DTest(const Array<int>& nx, const BasisFamily& basis) 
00138     : LP1DTestBase(nx), basis_(basis) {}
00139 
00140   /** */
00141   std::string name() const {
00142     TeuchosOStringStream ss;
00143     ss << "Poisson1D(basis=";
00144     basis_.print(ss);
00145     ss << ")";
00146     return ss.str();
00147   }
00148 
00149   /** */
00150   Array<int> pExpected() const {return tuple<int>(basis_.order()+1);}
00151 
00152   /** */
00153   Expr exactSoln() const 
00154     {
00155       Expr x = coord(0);
00156       const double pi = 4.0*atan(1.0);
00157 
00158       return sin(3.0*pi*x/2.0);
00159     }
00160 
00161   /** */
00162   LinearProblem prob(const Mesh& mesh) const
00163     {
00164       CellFilter left = domain().left();
00165 
00166       Expr u = new UnknownFunction(basis_, "u");
00167       Expr v = new TestFunction(basis_, "v");
00168       Expr dx = gradient(1);
00169       Expr x = coord(0);
00170 
00171       QuadratureFamily quad = new GaussianQuadrature(2*basis_.order());
00172 
00173       const double pi = 4.0*atan(1.0);
00174       Expr eqn = Integral(interior(), (dx*v)*(dx*u), quad)
00175         + Integral(interior(), -2.25*pi*pi*sin(3.0*pi*x/2.0)*v, quad);
00176 
00177       Expr bc = EssentialBC(left, v*u, quad);
00178       
00179       
00180       return LinearProblem(mesh, eqn, bc, v, u, vecType());
00181     }
00182 
00183   /** */
00184   Array<LPTestSpec> specs() const
00185     {
00186       double tol = 0.05;
00187       return tuple(
00188         LPTestSpec("amesos.xml", tol, makeSet<int>(1)),
00189         LPTestSpec("aztec-ifpack.xml", tol),
00190         LPTestSpec("aztec-ml.xml", tol),
00191         LPTestSpec("belos-ml.xml", tol),
00192         LPTestSpec("bicgstab.xml", tol)
00193         );
00194     }
00195 private:
00196   BasisFamily basis_;
00197 };
00198 
00199 
00200 
00201 
00202 /** Performs L2 projection in 1D with user-specified basis functions */
00203 class Projection1DTest : public LP1DTestBase
00204 {
00205 public:
00206   /** */
00207   Projection1DTest(const Array<int>& nx, const BasisFamily& basis) 
00208     : LP1DTestBase(nx), basis_(basis) {}
00209 
00210   /** */
00211   std::string name() const 
00212     {
00213       TeuchosOStringStream ss;
00214       ss << "Projection1D(basis=";
00215       basis_.print(ss);
00216       ss << ")";
00217       return ss.str();
00218     }
00219 
00220   /** */
00221   Array<int> pExpected() const {return tuple<int>(basis_.order()+1);}
00222 
00223   /** */
00224   Expr exactSoln() const 
00225     {
00226       Expr x = coord(0);
00227       return 10.0*x*exp(-1.0*x);
00228     }
00229 
00230   /** */
00231   LinearProblem prob(const Mesh& mesh) const
00232     {
00233       Expr u = new UnknownFunction(basis_, "u");
00234       Expr v = new TestFunction(basis_, "v");
00235       Expr x = coord(0);
00236 
00237       int p = basis_.order();
00238       QuadratureFamily quad = new GaussianQuadrature(2*p);
00239 
00240       Expr ex = exactSoln();
00241       Expr eqn = Integral(interior(), v*(u-ex), quad);
00242       Expr bc;
00243       
00244       return LinearProblem(mesh, eqn, bc, v, u, vecType());
00245     }
00246 
00247   /** */
00248   Array<LPTestSpec> specs() const
00249     {
00250       double tol = 0.05;
00251       return tuple(
00252         LPTestSpec("amesos.xml", tol, makeSet<int>(1)),
00253         LPTestSpec("aztec-ifpack.xml", tol),
00254         LPTestSpec("aztec-ml.xml", tol),
00255         LPTestSpec("bicgstab.xml", tol)
00256         );
00257     }
00258 private:
00259   BasisFamily basis_;
00260 };
00261 
00262 
00263 
00264 /** Solves the Helmholtz equation in 1D with Lagrange(2) basis functions */
00265 class Helmholtz1DTest : public LP1DTestBase
00266 {
00267 public:
00268   /** */
00269   Helmholtz1DTest(const Array<int>& nx) 
00270     : LP1DTestBase(0.0, atan(1.0), nx) {}
00271 
00272   /** */
00273   std::string name() const {return "Helmholtz1D";}
00274 
00275   /** */
00276   Expr exactSoln() const 
00277     {
00278       Expr x = coord(0);
00279       return cos(x) + sin(x);
00280     }
00281 
00282   /** */
00283   Array<int> pExpected() const {return tuple(3);}
00284 
00285   /** */
00286   LinearProblem prob(const Mesh& mesh) const
00287     {
00288       CellFilter left = domain().left();
00289       CellFilter right = domain().right();
00290 
00291       Expr u = new UnknownFunction(new Lagrange(2), "u");
00292       Expr v = new TestFunction(new Lagrange(2), "v");
00293       Expr dx = gradient(1);
00294       Expr x = coord(0);
00295 
00296       QuadratureFamily quad = new GaussianQuadrature(4);
00297 
00298       Expr eqn = Integral(interior(), (dx*v)*(dx*u) - v*u, quad);
00299       Expr bc = EssentialBC(left, v*(u-cos(x)), quad);
00300       
00301       return LinearProblem(mesh, eqn, bc, v, u, vecType());
00302     }
00303 
00304   /** */
00305   Array<LPTestSpec> specs() const
00306     {
00307       double tol = 0.05;
00308       return tuple(
00309         LPTestSpec("amesos.xml", tol, makeSet<int>(1)),
00310         LPTestSpec("aztec-ifpack.xml", tol),
00311         LPTestSpec("aztec-ml.xml", tol),
00312         LPTestSpec("belos-ml.xml", tol),
00313         LPTestSpec("bicgstab.xml", tol)
00314         );
00315     }
00316 private:
00317 };
00318 
00319 
00320 /** */
00321 class Coupled1DTest : public LP1DTestBase
00322 {
00323 public:
00324   /** */
00325   Coupled1DTest(const Array<int>& nx) : LP1DTestBase(nx) {}
00326 
00327   /** */
00328   std::string name() const {return "Coupled1D";}
00329 
00330   /** */
00331   Expr exactSoln() const 
00332     {
00333       Expr x = coord(0);
00334 
00335       Expr x2 = x*x;
00336       Expr x3 = x*x2;
00337       
00338       Expr uExact = (1.0/120.0)*x2*x3 - 1.0/36.0 * x3 + 7.0/360.0 * x;
00339       Expr vExact = 1.0/6.0 * x * (x2 - 1.0);
00340 
00341       return List(vExact, uExact);
00342     }
00343 
00344   /** */
00345   Array<int> pExpected() const {return tuple(2, 2);}
00346 
00347   /** */
00348   LinearProblem prob(const Mesh& mesh) const
00349     {
00350       CellFilter left = domain().left();
00351       CellFilter right = domain().right();
00352 
00353 
00354       Expr u = new UnknownFunction(new Lagrange(3), "u");
00355       Expr v = new UnknownFunction(new Lagrange(1), "v");
00356       Expr du = new TestFunction(new Lagrange(3), "du");
00357       Expr dv = new TestFunction(new Lagrange(1), "dv");
00358 
00359       Expr dx = gradient(1);
00360       Expr x = coord(0);
00361 
00362 
00363       QuadratureFamily quad = new GaussianQuadrature(10);
00364 
00365       Expr eqn = Integral(interior(), 
00366         (dx*du)*(dx*u) + du*v + (dx*dv)*(dx*v) + x*dv, 
00367         quad);
00368 
00369       Expr bc = EssentialBC(left, du*u + dv*v, quad)
00370         + EssentialBC(right, du*u + dv*v, quad);
00371 
00372 
00373       return LinearProblem(mesh, eqn, bc, 
00374         List(dv,du), List(v,u), vecType());
00375       
00376     }  
00377 
00378   /** */
00379   Array<LPTestSpec> specs() const
00380     {
00381       double tol = 0.05;
00382       return tuple(
00383         LPTestSpec("amesos.xml", tol, makeSet<int>(1)),
00384         LPTestSpec("aztec-ifpack.xml", tol),
00385         LPTestSpec("belos-ifpack.xml", tol),
00386         LPTestSpec("bicgstab.xml", tol)
00387         );
00388     }
00389 };
00390 
00391 
00392 
00393 
00394 
00395 
00396 
00397 int main(int argc, char** argv)
00398 {
00399   try
00400   {
00401     Sundance::init(&argc, &argv);
00402     Tabs::showDepth() = false;
00403     LinearSolveDriver::solveFailureIsFatal() = false;
00404     int np = MPIComm::world().getNProc();
00405 
00406     LPTestSuite tests;
00407 
00408     Array<int> nx = tuple(8, 16, 24, 32, 40);
00409 
00410     tests.registerTest(rcp(new SimplePoisson1DTest(nx)));
00411 
00412     for (int p=1; p<=3; p++)
00413     {
00414       if (np==1 || np % 4 == 0) 
00415       {
00416         tests.registerTest(rcp(new Poisson1DTest(nx, new Lagrange(p))));
00417         tests.registerTest(rcp(new Poisson1DTest(nx, new Bernstein(p))));
00418         tests.registerTest(rcp(new Projection1DTest(nx, new Lagrange(p))));
00419         tests.registerTest(rcp(new Projection1DTest(nx, new Bernstein(p))));
00420       }
00421       else
00422       {
00423         Out::root() << "skipping 1D tests needing numprocs=1 or 4" << std::endl;
00424       }
00425     }
00426 
00427     tests.registerTest(rcp(new Helmholtz1DTest(nx))); 
00428 
00429     
00430     tests.registerTest(rcp(new Coupled1DTest(nx))); 
00431 
00432     bool pass = tests.run();
00433 
00434     Out::root() << "total test status: " << pass << std::endl;
00435 
00436     Sundance::passFailTest(pass);
00437   }
00438   catch(std::exception& e)
00439   {
00440     Sundance::handleException(e);
00441   }
00442   Sundance::finalize(); return Sundance::testStatus(); 
00443 }

Site Contact