VorticityStokes2D.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 "Sundance.hpp"
00032 #include "SundanceEvaluator.hpp"
00033 #include "Teuchos_XMLParameterListWriter.hpp"
00034 
00035 using Sundance::List;
00036 /** 
00037  * Solves the 2D Stokes driven cavity in the vorticity-streamfunction 
00038  * formulation
00039  */
00040 
00041 CELL_PREDICATE(LeftPointTest, {return fabs(x[0]) < 1.0e-10;})
00042 CELL_PREDICATE(BottomPointTest, {return fabs(x[1]) < 1.0e-10;})
00043 CELL_PREDICATE(RightPointTest, {return fabs(x[0]-1.0) < 1.0e-10;})
00044 CELL_PREDICATE(TopPointTest, {return fabs(x[1]-1.0) < 1.0e-10;})
00045 
00046 
00047 int main(int argc, char** argv)
00048 {
00049   
00050   try
00051     {
00052       Sundance::init(&argc, &argv);
00053       int np = MPIComm::world().getNProc();
00054 
00055       /* We will do our linear algebra using Epetra */
00056       VectorType<double> vecType = new EpetraVectorType();
00057 
00058       /* Create a mesh. It will be of type BasisSimplicialMesh, and will
00059        * be built using a PartitionedRectangleMesher. */
00060       int nx = 32;
00061       int ny = 32;
00062       MeshType meshType = new BasicSimplicialMeshType();
00063       MeshSource mesher = new PartitionedRectangleMesher(0.0, 1.0, nx, np,
00064                                                          0.0, 1.0, ny, 1,
00065                                                          meshType);
00066       Mesh mesh = mesher.getMesh();
00067 
00068       /* Create a cell filter that will identify the maximal cells
00069        * in the interior of the domain */
00070       CellFilter interior = new MaximalCellFilter();
00071       CellFilter edges = new DimensionalCellFilter(1);
00072 
00073       CellFilter left = edges.coordSubset(0,0.0);
00074       CellFilter right = edges.coordSubset(0,1.0);
00075       CellFilter top = edges.subset(new TopPointTest());
00076       CellFilter bottom = edges.subset(new BottomPointTest());
00077 
00078       
00079       /* Create unknown and test functions, discretized using first-order
00080        * Lagrange interpolants */
00081       Expr psi = new UnknownFunction(new Lagrange(1), "psi");
00082       Expr vPsi = new TestFunction(new Lagrange(1), "vPsi");
00083       Expr omega = new UnknownFunction(new Lagrange(1), "omega");
00084       Expr vOmega = new TestFunction(new Lagrange(1), "vOmega");
00085 
00086       /* Create differential operator and coordinate functions */
00087       Expr dx = new Derivative(0);
00088       Expr dy = new Derivative(1);
00089       Expr grad = List(dx, dy);
00090       Expr x = new CoordExpr(0);
00091       Expr y = new CoordExpr(1);
00092 
00093       /* We need a quadrature rule for doing the integrations */
00094       QuadratureFamily quad2 = new GaussianQuadrature(2);
00095       QuadratureFamily quad4 = new GaussianQuadrature(4);
00096 
00097       /* Define the weak form */
00098       Expr eqn = Integral(interior, (grad*vPsi)*(grad*psi) 
00099                           + (grad*vOmega)*(grad*omega) + vPsi*omega, quad2)
00100         + Integral(top, -1.0*vPsi, quad4);
00101       /* Define the Dirichlet BC */
00102       Expr bc = EssentialBC(bottom, vOmega*psi, quad2) 
00103         + EssentialBC(top, vOmega*psi, quad2) 
00104         + EssentialBC(left, vOmega*psi, quad2) 
00105         + EssentialBC(right, vOmega*psi, quad2);
00106 
00107 
00108 
00109       /* We can now set up the linear problem! */
00110       LinearProblem prob(mesh, eqn, bc, List(vPsi, vOmega), 
00111                          List(psi, omega), vecType);
00112 
00113 
00114       LinearSolver<double> solver 
00115         = LinearSolverBuilder::createSolver("bicgstab.xml");
00116 
00117 
00118       std::cerr << "starting solve..." << std::endl;
00119 
00120       Expr soln = prob.solve(solver);
00121 
00122       /* Write the field in VTK format */
00123       FieldWriter w = new VTKWriter("VorticityStokes2D");
00124       w.addMesh(mesh);
00125       w.addField("psi", new ExprFieldWrapper(soln[0]));
00126       w.addField("omega", new ExprFieldWrapper(soln[1]));
00127       w.write();
00128 
00129       /* As a check, we integrate the vorticity over the domain. By 
00130        * Stokes' theorem this should be equal to the line integral
00131        * of the velocity around the boundary (which is 1). */
00132       Expr totalVorticityExpr = Integral(interior, soln[1], quad2);
00133       double totalVorticity = evaluateIntegral(mesh, totalVorticityExpr);
00134       std::cerr << "total vorticity = " << totalVorticity << std::endl;
00135 
00136       double tol = 1.0e-4;
00137       Sundance::passFailTest(fabs(totalVorticity-1.0), tol);
00138     }
00139   catch(std::exception& e)
00140     {
00141       std::cerr << e.what() << std::endl;
00142     }
00143   Sundance::finalize(); return Sundance::testStatus(); 
00144 }

Site Contact