VorticityStokes2D.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 #include "Teuchos_XMLParameterListWriter.hpp"
00045 
00046 using Sundance::List;
00047 /** 
00048  * Solves the 2D Stokes driven cavity in the vorticity-streamfunction 
00049  * formulation
00050  */
00051 
00052 CELL_PREDICATE(LeftPointTest, {return fabs(x[0]) < 1.0e-10;})
00053 CELL_PREDICATE(BottomPointTest, {return fabs(x[1]) < 1.0e-10;})
00054 CELL_PREDICATE(RightPointTest, {return fabs(x[0]-1.0) < 1.0e-10;})
00055 CELL_PREDICATE(TopPointTest, {return fabs(x[1]-1.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 = 32;
00072       int ny = 32;
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 edges = new DimensionalCellFilter(1);
00083 
00084       CellFilter left = edges.coordSubset(0,0.0);
00085       CellFilter right = edges.coordSubset(0,1.0);
00086       CellFilter top = edges.subset(new TopPointTest());
00087       CellFilter bottom = edges.subset(new BottomPointTest());
00088 
00089       
00090       /* Create unknown and test functions, discretized using first-order
00091        * Lagrange interpolants */
00092       Expr psi = new UnknownFunction(new Lagrange(1), "psi");
00093       Expr vPsi = new TestFunction(new Lagrange(1), "vPsi");
00094       Expr omega = new UnknownFunction(new Lagrange(1), "omega");
00095       Expr vOmega = new TestFunction(new Lagrange(1), "vOmega");
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       Expr eqn = Integral(interior, (grad*vPsi)*(grad*psi) 
00110                           + (grad*vOmega)*(grad*omega) + vPsi*omega, quad2)
00111         + Integral(top, -1.0*vPsi, quad4);
00112       /* Define the Dirichlet BC */
00113       Expr bc = EssentialBC(bottom, vOmega*psi, quad2) 
00114         + EssentialBC(top, vOmega*psi, quad2) 
00115         + EssentialBC(left, vOmega*psi, quad2) 
00116         + EssentialBC(right, vOmega*psi, quad2);
00117 
00118 
00119 
00120       /* We can now set up the linear problem! */
00121       LinearProblem prob(mesh, eqn, bc, List(vPsi, vOmega), 
00122                          List(psi, omega), vecType);
00123 
00124 
00125       LinearSolver<double> solver 
00126         = LinearSolverBuilder::createSolver("bicgstab.xml");
00127 
00128 
00129       std::cerr << "starting solve..." << std::endl;
00130 
00131       Expr soln = prob.solve(solver);
00132 
00133       /* Write the field in VTK format */
00134       FieldWriter w = new VTKWriter("VorticityStokes2D");
00135       w.addMesh(mesh);
00136       w.addField("psi", new ExprFieldWrapper(soln[0]));
00137       w.addField("omega", new ExprFieldWrapper(soln[1]));
00138       w.write();
00139 
00140       /* As a check, we integrate the vorticity over the domain. By 
00141        * Stokes' theorem this should be equal to the line integral
00142        * of the velocity around the boundary (which is 1). */
00143       Expr totalVorticityExpr = Integral(interior, soln[1], quad2);
00144       double totalVorticity = evaluateIntegral(mesh, totalVorticityExpr);
00145       std::cerr << "total vorticity = " << totalVorticity << std::endl;
00146 
00147       double tol = 1.0e-4;
00148       Sundance::passFailTest(fabs(totalVorticity-1.0), tol);
00149     }
00150   catch(std::exception& e)
00151     {
00152       std::cerr << e.what() << std::endl;
00153     }
00154   Sundance::finalize(); return Sundance::testStatus(); 
00155 }

Site Contact