TransientNonlinear1D.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 initial-boundary-value problem in 
00047  * one spatial dimension using Galerkin finite elements in space and
00048  * Crank-Nicolson stepping in time. 
00049  *
00050  * ------------------------------------------------------------------------ */
00051 
00052 
00053 
00054 Expr force(const double& epsilon, const Expr& x, const Expr& t)
00055 {
00056   const double pi = 4.0*atan(1.0);
00057 
00058   Expr rtn = -8.0*epsilon*cos(2.0*pi*t)*
00059     pow(1.0 + pow(x,2.0)*epsilon*cos(2.0*pi*t),2.0)*
00060     (1.0 + 7.0*pow(x,2.0)*epsilon*cos(2.0*pi*t)) - 
00061     2.0*pi*pow(x,2.0)*epsilon*sin(2.0*pi*t);
00062 
00063   return rtn;
00064 }
00065 
00066 
00067 int main(int argc, char** argv)
00068 {
00069   const double pi = 4.0*atan(1.0);
00070 
00071   try
00072   {
00073     Sundance::init(&argc, &argv);
00074     int np = MPIComm::world().getNProc();
00075 
00076     /* We will do our linear algebra using Epetra */
00077     VectorType<double> vecType = new EpetraVectorType();
00078 
00079     /* Create a mesh. It will be of type BasisSimplicialMesh, and will
00080      * be built using a PartitionedLineMesher. */
00081     MeshType meshType = new BasicSimplicialMeshType();
00082     int nx = 32;
00083     MeshSource mesher = new PartitionedLineMesher(0.0, 1.0, nx*np, meshType);
00084     Mesh mesh = mesher.getMesh();
00085 
00086     /* Create a cell filter that will identify the maximal cells
00087      * in the interior of the domain */
00088     CellFilter interior = new MaximalCellFilter();
00089     CellFilter points = new DimensionalCellFilter(0);
00090     CellFilter leftPoint = points.coordSubset(0, 0.0);
00091     CellFilter rightPoint = points.coordSubset(0, 1.0);
00092 
00093       
00094     /* Create unknown and test functions, discretized using first-order
00095      * Lagrange interpolants */
00096     BasisFamily bas = new Lagrange(1);
00097     Expr u = new UnknownFunction(bas, "u");
00098     Expr v = new TestFunction(bas, "v");
00099 
00100     /* Create differential operator and coordinate function */
00101     Expr dx = new Derivative(0);
00102     Expr x = new CoordExpr(0);
00103 
00104     double epsilon = 0.25;
00105     /* The initial profile is u(x,0)=1 + epsilon x^2. 
00106      * Project this onto a discrete function. This discrete function,
00107      * uPrev, will represent
00108      * the initial conditions but will also be re-used as the starting value for
00109      * each timestep */
00110     Expr uStart = 1.0 + epsilon*x*x;
00111     DiscreteSpace discSpace(mesh, bas, vecType);
00112     L2Projector projector(discSpace, 1.0 + epsilon*x*x);
00113     Expr uPrev = projector.project();
00114 
00115     /* We need another discrete function for the current Newton approximation */
00116     Expr uNewt = copyDiscreteFunction(uPrev, "uNewt");
00117 
00118     /* We need a quadrature rule for doing the integrations */
00119     QuadratureFamily quad = new GaussianQuadrature(4);
00120 
00121     int nSteps = nx;
00122     double dt = 1.0/((double) nSteps);
00123     /* Represent the time variable as a parameter expression, NOT as
00124      * a double variable. The reason is that we need to be able to update
00125      * the time value without rebuilding expressions. */
00126     Expr t = new Sundance::Parameter(0.0);
00127     Expr tPrev = new Sundance::Parameter(0.0);
00128 
00129     /* Define the weak form, semidiscretized in time */
00130     Expr eqn = Integral(interior, v*(u-uPrev) 
00131       + dt/2.0*(dx*v)*((dx*pow(u, 4.0))+(dx*pow(uPrev, 4.0)))
00132       - dt/2.0*v*(force(epsilon, x, t)+force(epsilon, x, tPrev)), quad); 
00133     /* Define the Dirichlet BC on the right side of the domain */
00134     Expr bc = EssentialBC(rightPoint, v*(u - 1.0 - epsilon*cos(2.0*pi*t)),quad);
00135 
00136     /* We can now set up the nonlinear problem! */
00137     NonlinearProblem prob(mesh, eqn, bc, v, u, uNewt, vecType); 
00138 
00139     NonlinearSolver<double> solver 
00140       = NonlinearSolverBuilder::createSolver("playa-newton-amesos.xml");
00141 
00142     /* Wrtie the initial conditions */
00143     FieldWriter w = new DSVWriter("transientNonlinear1D-0.dat"); 
00144     w.addMesh(mesh);
00145     w.addField("u", new ExprFieldWrapper(uPrev[0]));
00146     w.write();
00147 
00148     /* loop over timesteps */
00149     for (int i=0; i<nSteps; i++)
00150     {
00151       /* Set the times t_i and t_{i+1} */
00152       Out::root() << "timestep #" << i << endl;
00153       t.setParameterValue((i+1)*dt);
00154       tPrev.setParameterValue(i*dt);
00155 
00156       SolverState<double> state = prob.solve(solver);
00157 
00158       TEUCHOS_TEST_FOR_EXCEPTION(state.finalState() != SolveConverged,
00159         std::runtime_error,
00160         "Nonlinear solve failed to converge: message=" << state.finalMsg());
00161           
00162       updateDiscreteFunction(uNewt, uPrev);
00163       
00164       FieldWriter writer = new MatlabWriter("transientNonlinear1D-" 
00165         + Teuchos::toString(i+1) + ".dat");
00166       writer.addMesh(mesh);
00167       writer.addField("u", new ExprFieldWrapper(uPrev[0]));
00168       writer.write();
00169     }
00170     
00171     double err = L2Norm(mesh, interior, uPrev - uStart, quad);
00172     Out::root() << "dt=" << setw(16) << dt << " error = " << setw(16) << err << endl;
00173     
00174     double h = 1.0/(nx-1.0);
00175     double tol = 0.1*(dt*dt + h*h);
00176     Sundance::passFailTest(err, tol);
00177   }
00178   catch(std::exception& e)
00179   {
00180     Sundance::handleException(e);
00181   }
00182   Sundance::finalize(); return Sundance::testStatus(); 
00183 }

Site Contact