RadDiff.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 
00044 /* ------------------------------------------------------------------------ 
00045  *
00046  * This program solves a nonlinear boundary-value problem in two spatial
00047  * dimensions using Galerkin finite elements and Newton's method
00048  *
00049  * ------------------------------------------------------------------------ */
00050 
00051 const double pi = 4.0*atan(1.0);
00052 
00053 
00054 
00055 class WestEdgeTest : public CellPredicateFunctorBase, 
00056                      public Playa::Handleable<CellPredicateFunctorBase>
00057 {
00058 public:
00059   WestEdgeTest() : CellPredicateFunctorBase("WestEdgeTest") {}
00060   virtual ~WestEdgeTest() {}
00061   virtual bool operator()(const Point& x) const {return fabs(x[0])<1.0e-10;}
00062   GET_RCP(CellPredicateFunctorBase);
00063 };
00064 
00065 
00066 CELL_PREDICATE(EastEdgeTest, {return fabs(x[0]-1.0) < 1.0e-10;})
00067 CELL_PREDICATE(SouthEdgeTest, {return fabs(x[1]) < 1.0e-10;})
00068 CELL_PREDICATE(NorthEdgeTest, {return fabs(x[1]-1.0) < 1.0e-10;})
00069 
00070 bool inlineNewtonSolve(NonlinearProblem prob,
00071   Expr uNewt,
00072   int maxIters, 
00073   double newtTol);
00074 
00075 bool noxNewtonSolve(const NonlinearProblem& prob, Expr uNewt);
00076 
00077 bool playaNewtonArmijoSolve(const NonlinearProblem& prob, Expr uNewt);
00078 
00079 int main(int argc, char** argv)
00080 {
00081   try
00082     {
00083       Sundance::init(&argc, &argv);
00084 
00085       /* We will do our linear algebra using Epetra */
00086       VectorType<double> vecType = new EpetraVectorType();
00087 
00088       /* Create a mesh. It will be of type BasisSimplicialMesh, and will
00089        * be built using a PartitionedLineMesher. */
00090       MeshType meshType = new BasicSimplicialMeshType();
00091       int nx = 128;
00092       MeshSource mesher = new PartitionedRectangleMesher(
00093         0.0, 1.0, nx, 1, 
00094         0.0, 1.0, nx, 1, meshType);
00095       Mesh mesh = mesher.getMesh();
00096 
00097       /* Create a cell filter that will identify the maximal cells
00098        * in the interior of the domain */
00099       CellFilter interior = new MaximalCellFilter();
00100       CellFilter edges = new DimensionalCellFilter(1);
00101       CellFilter northEdge = edges.subset(new NorthEdgeTest());
00102       CellFilter southEdge = edges.subset(new SouthEdgeTest());
00103       CellFilter eastEdge = edges.subset(new EastEdgeTest());
00104       CellFilter westEdge = edges.subset(new WestEdgeTest());
00105 
00106       
00107       /* Create unknown and test functions, discretized using first-order
00108        * Lagrange interpolants */
00109       BasisFamily bas = new Lagrange(1);
00110       Expr u = new UnknownFunction(bas, "u");
00111       Expr v = new TestFunction(bas, "v");
00112 
00113       /* Create differential operator and coordinate function */
00114       Expr dx = new Derivative(0);
00115       Expr dy = new Derivative(1);
00116       Expr grad = List(dx, dy);
00117       Expr x = new CoordExpr(0);
00118       Expr y = new CoordExpr(1);
00119 
00120       /* The initial guess is u(x,y)=1 */
00121       DiscreteSpace discSpace(mesh, bas, vecType);
00122       Expr uNewt = new DiscreteFunction(discSpace, 1.0, "uNewt");
00123 
00124       /* We need a quadrature rule for doing the integrations */
00125       QuadratureFamily quad = new GaussianQuadrature(4);
00126  WatchFlag watchMe("watch me");
00127     watchMe.setParam("symbolic preprocessing", 1);
00128     watchMe.setParam("discrete function evaluation", 3);
00129     watchMe.setParam("integration setup", 6);
00130     watchMe.setParam("integration", 6);
00131     watchMe.setParam("fill", 6);
00132     watchMe.setParam("evaluation", 2);
00133     watchMe.deactivate();  // Deactive
00134       /* Define the weak form */
00135       Expr eqn = Integral(interior, u*u*u*(grad*v)*(grad*u), quad);
00136       /* Define the Dirichlet BC */
00137       Expr bc = EssentialBC(northEdge, v*(u - pow(1.0 + sin(pi*x),0.25)),quad)
00138         + EssentialBC(southEdge + eastEdge + westEdge, 
00139           v*(u - 1.0),quad);
00140 
00141       /* We can now set up the nonlinear problem! */
00142       NonlinearProblem prob(mesh, eqn, bc, v, u, uNewt, vecType); 
00143 
00144       bool nonlinSolveOK = false;
00145       bool useNox = false;
00146       if (useNox)
00147       {
00148         nonlinSolveOK = noxNewtonSolve(prob, uNewt);
00149       }
00150       else
00151       {
00152         nonlinSolveOK = playaNewtonArmijoSolve(prob, uNewt);
00153       }
00154       
00155       TEUCHOS_TEST_FOR_EXCEPT(!nonlinSolveOK);
00156 
00157       FieldWriter writer = new VTKWriter("SteadyRadDiff2D"); 
00158       writer.addMesh(mesh);
00159       writer.addField("u", new ExprFieldWrapper(uNewt[0]));
00160       writer.write();
00161 
00162       Expr uExact = pow(1.0 + sin(pi*x)*sinh(pi*y)/sinh(pi), 0.25);
00163       double err = L2Norm(mesh, interior, uNewt-uExact, quad);
00164       Out::root() << "error = " << setw(16) << err << endl;
00165 
00166       double tol = 1.0e-4;
00167       Sundance::passFailTest(err, tol);
00168     }
00169   catch(std::exception& e)
00170     {
00171       Sundance::handleException(e);
00172     }
00173   Sundance::finalize(); return Sundance::testStatus(); 
00174 }
00175 
00176 
00177 bool inlineNewtonSolve(NonlinearProblem prob,
00178   Expr uNewt,
00179   int maxIters, 
00180   double newtTol)
00181 {
00182   bool newtonConverged = false;
00183   LinearSolver<double> linSolver 
00184     = LinearSolverBuilder::createSolver("amesos.xml");
00185 
00186   /* Allocate objects for the Jacobian, residual, and Newton step */
00187   LinearOperator<double> J = prob.allocateJacobian();
00188   Vector<double> resid = J.range().createMember();
00189   Vector<double> newtonStep = J.domain().createMember();
00190 
00191   for (int i=0; i<maxIters; i++)
00192   {
00193     prob.setEvalPoint(uNewt);
00194     prob.computeJacobianAndFunction(J, resid);
00195     SolverState<double> solveState 
00196       = linSolver.solve(J, -1.0*resid, newtonStep);
00197     
00198     TEUCHOS_TEST_FOR_EXCEPTION(solveState.finalState() != SolveConverged,
00199       std::runtime_error,
00200       "linear solve failed!");
00201     
00202     addVecToDiscreteFunction(uNewt, newtonStep);
00203     double newtStepNorm = newtonStep.norm2();
00204     Out::root() << "|newt step| = " << newtStepNorm << endl;
00205     if (newtStepNorm < newtTol) 
00206     {
00207       newtonConverged = true;
00208       break;
00209     }
00210   }
00211 
00212   return newtonConverged;
00213 }
00214 
00215 bool noxNewtonSolve(const NonlinearProblem& prob, Expr uNewt)
00216 {
00217   /* Use a NOX nonlinear solver */
00218   ParameterXMLFileReader reader("nox-amesos.xml");
00219   ParameterList solverParams = reader.getParameters();
00220   NOXSolver solver(solverParams);
00221   
00222   SolverState<double> stat = prob.solve(solver);
00223 
00224   return stat.finalState()==SolveConverged;
00225 }
00226 
00227 bool playaNewtonArmijoSolve(const NonlinearProblem& prob, Expr uNewt)
00228 {
00229   /* Use the Playa Newton-Armijo nonlinear solver */
00230   LinearSolver<double> linSolver = LinearSolverBuilder::createSolver("amesos.xml");
00231   ParameterList solverParams("NewtonArmijoSolver");
00232   solverParams.set("Tau Relative", 1.0e-12);
00233   solverParams.set("Tau Absolute", 1.0e-12);
00234   solverParams.set("Alpha", 1.0e-4);
00235   solverParams.set("Verbosity", 3);
00236 
00237   NonlinearSolver<double> solver = new NewtonArmijoSolver<double>(solverParams, linSolver);
00238   
00239   SolverState<double> stat = prob.solve(solver);
00240   if (stat.finalState() != SolveConverged)
00241   {
00242     Out::os() << stat.finalMsg() << endl;
00243   }
00244 
00245   return stat.finalState()==SolveConverged;
00246 }

Site Contact