TaylorHoodStokes2D.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 "SundanceEvaluator.hpp"
00044 
00045 using Sundance::List;
00046 /** 
00047  * Solves the Poisson equation in 2D
00048  */
00049 
00050 CELL_PREDICATE(LeftPointTest, {return fabs(x[0]) < 1.0e-10;})
00051 CELL_PREDICATE(BottomPointTest, {return fabs(x[1]) < 1.0e-10;})
00052 CELL_PREDICATE(RightPointTest, {return fabs(x[0]-1.0) < 1.0e-10;})
00053 CELL_PREDICATE(TopPointTest, {return fabs(x[1]-1.0) < 1.0e-10;})
00054 
00055 CELL_PREDICATE(CornerPointTest, {return fabs(x[1]) < 1.0e-10 && fabs(x[0])<1.0e-10;})
00056 
00057 
00058 int main(int argc, char** argv)
00059 {
00060   
00061   try
00062     {
00063       Sundance::init(&argc, &argv);
00064       int np = MPIComm::world().getNProc();
00065 
00066       /* We will do our linear algebra using Epetra */
00067       VectorType<double> vecType = new EpetraVectorType();
00068 
00069       /* Create a mesh. It will be of type BasisSimplicialMesh, and will
00070        * be built using a PartitionedRectangleMesher. */
00071       int nx = 8;
00072       int ny = 8;
00073       MeshType meshType = new BasicSimplicialMeshType();
00074       MeshSource mesher = new PartitionedRectangleMesher(0.0, 1.0, nx, np,
00075                                                          0.0, 1.0, ny, 1,
00076                                                          meshType);
00077       Mesh mesh = mesher.getMesh();
00078 
00079       /* Create a cell filter that will identify the maximal cells
00080        * in the interior of the domain */
00081       CellFilter interior = new MaximalCellFilter();
00082       CellFilter bdry = new BoundaryCellFilter();
00083       CellFilter nodes = new DimensionalCellFilter(0);
00084       CellFilter corner = nodes.subset(new CornerPointTest());
00085 
00086       
00087       /* Unknown and test functions, using Taylor-Hood discretization */
00088       Expr ux = new UnknownFunction(new Lagrange(2), "u_x");
00089       Expr vx = new TestFunction(new Lagrange(2), "v_x");
00090       Expr uy = new UnknownFunction(new Lagrange(2), "u_y");
00091       Expr vy = new TestFunction(new Lagrange(2), "v_y");
00092       Expr p = new UnknownFunction(new Lagrange(1), "p");
00093       Expr q = new TestFunction(new Lagrange(1), "q");
00094       Expr u = List(ux, uy);
00095       Expr v = List(vx, vy);
00096 
00097       /* Create differential operator and coordinate functions */
00098       Expr dx = new Derivative(0);
00099       Expr dy = new Derivative(1);
00100       Expr grad = List(dx, dy);
00101       Expr x = new CoordExpr(0);
00102       Expr y = new CoordExpr(1);
00103 
00104       /* We need a quadrature rule for doing the integrations */
00105       QuadratureFamily quad2 = new GaussianQuadrature(2);
00106       QuadratureFamily quad4 = new GaussianQuadrature(4);
00107 
00108       /* Define the weak form */
00109 
00110       double pi = 4.0*atan(1.0);
00111       Expr sx = sin(pi*x);
00112       Expr cx = cos(pi*x);
00113       Expr sy = sin(pi*y);
00114       Expr cy = cos(pi*y);
00115       Expr psiExact = pow(pi, -3.0) * sx*sy;
00116       Expr uExact = pow(pi, -2.0)*List(-sx*cy, cx*sy);
00117       Expr fy = 4.0*cx*sy;
00118 
00119       Expr eqn = Integral(interior, (grad*vx)*(grad*ux)  
00120                           + (grad*vy)*(grad*uy)  - p*(dx*vx+dy*vy)
00121                           + q*(dx*ux+dy*uy) - vy*fy,
00122                           quad2);
00123         
00124       /* Define the Dirichlet BC */
00125       Expr bc = EssentialBC(bdry, v*(u-uExact), quad4)
00126         + EssentialBC(corner, q*p, quad4);
00127 
00128       /* We can now set up the linear problem! */
00129       LinearProblem prob(mesh, eqn, bc, List(vx, vy, q),
00130                          List(ux, uy, p), vecType);
00131 
00132 
00133 #ifdef HAVE_CONFIG_H
00134       ParameterXMLFileReader reader(searchForFile("SolverParameters/aztec-native.xml"));
00135 #else
00136       ParameterXMLFileReader reader("aztec-native.xml");
00137 #endif
00138       ParameterList solverParams = reader.getParameters();
00139       LinearSolver<double> solver 
00140         = LinearSolverBuilder::createSolver(solverParams);
00141 
00142       Expr soln = prob.solve(solver);
00143 
00144 
00145       /* Write the field in VTK format */
00146       FieldWriter w = new VTKWriter("TaylorHoodStokes2d");
00147       w.addMesh(mesh);
00148       w.addField("ux", new ExprFieldWrapper(soln[0]));
00149       w.addField("uy", new ExprFieldWrapper(soln[1]));
00150       w.addField("p", new ExprFieldWrapper(soln[2]));
00151       w.write();
00152 
00153       Expr err = List(soln[0], soln[1]) - uExact;
00154       Expr errExpr = Integral(interior, 
00155                               err*err,
00156                               quad4);
00157 
00158       FunctionalEvaluator errInt(mesh, errExpr);
00159 
00160       double errorSq = errInt.evaluate();
00161       std::cerr << "velocity error norm = " << sqrt(errorSq) << std::endl << std::endl;
00162 
00163       Sundance::passFailTest(sqrt(errorSq), 1.0e-3);
00164 
00165     }
00166   catch(std::exception& e)
00167     {
00168       Sundance::handleException(e);
00169     }
00170   Sundance::finalize(); return Sundance::testStatus(); 
00171 }

Site Contact