Go to the documentation of this file.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
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054 Expr force(const double& epsilon, const Expr& x, const Expr& t)
00055 {
00056 const double pi = 4.0*atan(1.0);
00057
00058 Expr rtn = -8.0*epsilon*cos(2.0*pi*t)*
00059 pow(1.0 + pow(x,2.0)*epsilon*cos(2.0*pi*t),2.0)*
00060 (1.0 + 7.0*pow(x,2.0)*epsilon*cos(2.0*pi*t)) -
00061 2.0*pi*pow(x,2.0)*epsilon*sin(2.0*pi*t);
00062
00063 return rtn;
00064 }
00065
00066
00067 int main(int argc, char** argv)
00068 {
00069 const double pi = 4.0*atan(1.0);
00070
00071 try
00072 {
00073 Sundance::init(&argc, &argv);
00074 int np = MPIComm::world().getNProc();
00075
00076
00077 VectorType<double> vecType = new EpetraVectorType();
00078
00079
00080
00081 MeshType meshType = new BasicSimplicialMeshType();
00082 int nx = 32;
00083 MeshSource mesher = new PartitionedLineMesher(0.0, 1.0, nx*np, meshType);
00084 Mesh mesh = mesher.getMesh();
00085
00086
00087
00088 CellFilter interior = new MaximalCellFilter();
00089 CellFilter points = new DimensionalCellFilter(0);
00090 CellFilter leftPoint = points.coordSubset(0, 0.0);
00091 CellFilter rightPoint = points.coordSubset(0, 1.0);
00092
00093
00094
00095
00096 BasisFamily bas = new Lagrange(1);
00097 Expr u = new UnknownFunction(bas, "u");
00098 Expr v = new TestFunction(bas, "v");
00099
00100
00101 Expr dx = new Derivative(0);
00102 Expr x = new CoordExpr(0);
00103
00104 double epsilon = 0.25;
00105
00106
00107
00108
00109
00110 Expr uStart = 1.0 + epsilon*x*x;
00111 DiscreteSpace discSpace(mesh, bas, vecType);
00112 L2Projector projector(discSpace, 1.0 + epsilon*x*x);
00113 Expr uPrev = projector.project();
00114
00115
00116 Expr uNewt = copyDiscreteFunction(uPrev, "uNewt");
00117
00118
00119 QuadratureFamily quad = new GaussianQuadrature(4);
00120
00121 int nSteps = nx;
00122 double dt = 1.0/((double) nSteps);
00123
00124
00125
00126 Expr t = new Sundance::Parameter(0.0);
00127 Expr tPrev = new Sundance::Parameter(0.0);
00128
00129
00130 Expr eqn = Integral(interior, v*(u-uPrev)
00131 + dt/2.0*(dx*v)*((dx*pow(u, 4.0))+(dx*pow(uPrev, 4.0)))
00132 - dt/2.0*v*(force(epsilon, x, t)+force(epsilon, x, tPrev)), quad);
00133
00134 Expr bc = EssentialBC(rightPoint, v*(u - 1.0 - epsilon*cos(2.0*pi*t)),quad);
00135
00136
00137 NonlinearProblem prob(mesh, eqn, bc, v, u, uNewt, vecType);
00138
00139 NonlinearSolver<double> solver
00140 = NonlinearSolverBuilder::createSolver("playa-newton-amesos.xml");
00141
00142
00143 FieldWriter w = new DSVWriter("transientNonlinear1D-0.dat");
00144 w.addMesh(mesh);
00145 w.addField("u", new ExprFieldWrapper(uPrev[0]));
00146 w.write();
00147
00148
00149 for (int i=0; i<nSteps; i++)
00150 {
00151
00152 Out::root() << "timestep #" << i << endl;
00153 t.setParameterValue((i+1)*dt);
00154 tPrev.setParameterValue(i*dt);
00155
00156 SolverState<double> state = prob.solve(solver);
00157
00158 TEUCHOS_TEST_FOR_EXCEPTION(state.finalState() != SolveConverged,
00159 std::runtime_error,
00160 "Nonlinear solve failed to converge: message=" << state.finalMsg());
00161
00162 updateDiscreteFunction(uNewt, uPrev);
00163
00164 FieldWriter writer = new MatlabWriter("transientNonlinear1D-"
00165 + Teuchos::toString(i+1) + ".dat");
00166 writer.addMesh(mesh);
00167 writer.addField("u", new ExprFieldWrapper(uPrev[0]));
00168 writer.write();
00169 }
00170
00171 double err = L2Norm(mesh, interior, uPrev - uStart, quad);
00172 Out::root() << "dt=" << setw(16) << dt << " error = " << setw(16) << err << endl;
00173
00174 double h = 1.0/(nx-1.0);
00175 double tol = 0.1*(dt*dt + h*h);
00176 Sundance::passFailTest(err, tol);
00177 }
00178 catch(std::exception& e)
00179 {
00180 Sundance::handleException(e);
00181 }
00182 Sundance::finalize(); return Sundance::testStatus();
00183 }