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 #include "Sundance.hpp"
00032
00033 using Sundance::List;
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052 int main(int argc, char** argv)
00053 {
00054
00055 try
00056 {
00057 Sundance::init(&argc, &argv);
00058
00059
00060 VectorType<double> vecType = new EpetraVectorType();
00061
00062
00063
00064 MeshType meshType = new BasicSimplicialMeshType();
00065 int nx = 32;
00066
00067 MeshSource mesher = new PartitionedRectangleMesher(0.0, 1.0, nx,
00068 0.0, 2.0, nx,
00069 meshType);
00070 Mesh mesh = mesher.getMesh();
00071
00072
00073
00074 CellFilter interior = new MaximalCellFilter();
00075
00076
00077 CellFilter edges = new DimensionalCellFilter(1);
00078
00079 CellFilter west = edges.coordSubset(0, 0.0);
00080 CellFilter east = edges.coordSubset(0, 1.0);
00081
00082
00083
00084
00085 Expr u1 = new UnknownFunction(new Lagrange(3), "u");
00086 Expr u2 = new UnknownFunction(new Lagrange(2), "v");
00087 Expr du1 = new TestFunction(new Lagrange(3), "du");
00088 Expr du2 = new TestFunction(new Lagrange(2), "dv");
00089
00090
00091 Expr dx = new Derivative(0);
00092 Expr x = new CoordExpr(0);
00093 Expr dy = new Derivative(1);
00094 Expr y = new CoordExpr(1);
00095 Expr grad = List(dx, dy);
00096
00097
00098 QuadratureFamily quad = new GaussianQuadrature(6);
00099
00100
00101
00102 Expr eqn = Integral(interior,
00103 (grad*du1)*(grad*u1) + du1*u2 + (grad*du2)*(grad*u2) + x*du2,
00104 quad);
00105
00106 Expr bc = EssentialBC(east + west, du1*u1 + du2*u2, quad);
00107
00108
00109
00110
00111 LinearProblem prob(mesh, eqn, bc, List(du1,du2), List(u1,u2), vecType);
00112
00113
00114
00115
00116 LinearSolver<double> solver = LinearSolverBuilder::createSolver("aztec-ml.xml");
00117
00118
00119
00120 Expr soln = prob.solve(solver);
00121
00122
00123 Expr x2 = x*x;
00124 Expr x3 = x*x2;
00125
00126 Expr u1Exact = (1.0/120.0)*x2*x3 - 1.0/36.0 * x3 + 7.0/360.0 * x;
00127 Expr u2Exact = 1.0/6.0 * x * (x2 - 1.0);
00128
00129 Expr u1Err = u1Exact - soln[0];
00130 Expr u2Err = u2Exact - soln[1];
00131
00132 double u1ErrorNorm = L2Norm(mesh, interior, u1Err, quad);
00133 std::cerr << "u1 error norm = " << u1ErrorNorm << std::endl << std::endl;
00134
00135 double u2ErrorNorm = L2Norm(mesh, interior, u2Err, quad);
00136 std::cerr << "u2 error norm = " << u2ErrorNorm << std::endl << std::endl;
00137
00138
00139
00140
00141 FieldWriter w = new VTKWriter("Coupled2D");
00142 w.addMesh(mesh);
00143 w.addField("v", new ExprFieldWrapper(soln[0]));
00144 w.addField("u", new ExprFieldWrapper(soln[1]));
00145 w.write();
00146
00147
00148
00149 double tol = 1.0e-5;
00150 Sundance::passFailTest(max(u1ErrorNorm, u2ErrorNorm), tol);
00151 }
00152 catch(std::exception& e)
00153 {
00154 Sundance::handleException(e);
00155 }
00156 Sundance::finalize();
00157 return Sundance::testStatus();
00158 }