00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042 #include "Sundance.hpp"
00043 #include "SundanceEvaluator.hpp"
00044
00045 using Sundance::List;
00046
00047
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
00067 VectorType<double> vecType = new EpetraVectorType();
00068
00069
00070
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
00080
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
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
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
00105 QuadratureFamily quad2 = new GaussianQuadrature(2);
00106 QuadratureFamily quad4 = new GaussianQuadrature(4);
00107
00108
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
00125 Expr bc = EssentialBC(bdry, v*(u-uExact), quad4)
00126 + EssentialBC(corner, q*p, quad4);
00127
00128
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
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 }