75 virtual ~FE_Evolution() { }
79 int main(
int argc,
char *argv[])
83 const char *mesh_file =
"../data/periodic-hexagon.mesh";
86 int ode_solver_type = 4;
87 double t_final = 10.0;
89 bool visualization =
true;
94 cout.precision(precision);
97 args.
AddOption(&mesh_file,
"-m",
"--mesh",
100 "Problem setup to use. See options in velocity_function().");
101 args.
AddOption(&ref_levels,
"-r",
"--refine",
102 "Number of times to refine the mesh uniformly.");
104 "Order (degree) of the finite elements.");
105 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
106 "ODE solver: 1 - Forward Euler, 2 - RK2 SSP, 3 - RK3 SSP,"
107 " 4 - RK4, 6 - RK6.");
108 args.
AddOption(&t_final,
"-tf",
"--t-final",
109 "Final time; start time is 0.");
110 args.
AddOption(&dt,
"-dt",
"--time-step",
112 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
113 "--no-visualization",
114 "Enable or disable GLVis visualization.");
115 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
116 "--no-visit-datafiles",
117 "Save data files for VisIt (visit.llnl.gov) visualization.");
118 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
119 "Visualize every n-th timestep.");
130 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
136 switch (ode_solver_type)
139 case 2: ode_solver =
new RK2Solver(1.0);
break;
141 case 4: ode_solver =
new RK4Solver;
break;
142 case 6: ode_solver =
new RK6Solver;
break;
144 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
152 for (
int lev = 0; lev < ref_levels; lev++)
167 cout <<
"Number of unknowns: " << fes.
GetVSize() << endl;
203 ofstream omesh(
"ex9.mesh");
204 omesh.precision(precision);
206 ofstream osol(
"ex9-init.gf");
207 osol.precision(precision);
223 char vishost[] =
"localhost";
225 sout.
open(vishost, visport);
228 cout <<
"Unable to connect to GLVis server at "
229 << vishost <<
':' << visport << endl;
230 visualization =
false;
231 cout <<
"GLVis visualization disabled.\n";
235 sout.precision(precision);
236 sout <<
"solution\n" << *mesh << u;
239 cout <<
"GLVis visualization paused."
240 <<
" Press space (in the GLVis window) to resume it.\n";
248 ode_solver->
Init(adv);
251 for (
int ti = 0;
true; )
253 if (t >= t_final - dt/2)
258 ode_solver->
Step(u, t, dt);
261 if (ti % vis_steps == 0)
263 cout <<
"time step: " << ti <<
", time: " << t << endl;
267 sout <<
"solution\n" << *mesh << u << flush;
282 ofstream osol(
"ex9-final.gf");
283 osol.precision(precision);
299 M_solver.SetPreconditioner(M_prec);
300 M_solver.SetOperator(M);
302 M_solver.iterative_mode =
false;
303 M_solver.SetRelTol(1e-9);
304 M_solver.SetAbsTol(0.0);
305 M_solver.SetMaxIter(100);
306 M_solver.SetPrintLevel(0);
325 for (
int i = 0; i <
dim; i++)
337 case 1: v(0) = 1.0;
break;
338 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
339 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
348 const double w = M_PI/2;
351 case 1: v(0) = 1.0;
break;
352 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
353 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
360 const double w = M_PI/2;
361 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
365 case 1: v(0) = 1.0;
break;
366 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
367 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
380 for (
int i = 0; i <
dim; i++)
393 return exp(-40.*pow(X(0)-0.5,2));
397 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
400 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
404 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
405 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
411 double x_ = X(0), y_ = X(1), rho, phi;
412 rho = hypot(x_, y_) ;
414 return pow(sin(M_PI*rho),2)*sin(3*phi);
418 const double f = M_PI;
419 return sin(f*X(0))*sin(f*X(1));
Conjugate gradient method.
Class for grid function - Vector with associated FE space.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Data type for scaled Jacobi-type smoother of sparse matrix.
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Base abstract class for time dependent operators: (x,t) -> f(x,t)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
virtual void Step(Vector &x, double &t, double &dt)=0
int Size() const
Returns the size of the vector.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
virtual void Save()
Save the collection and a VisIt root file.
virtual void Init(TimeDependentOperator &_f)
int main(int argc, char *argv[])
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
void PrintUsage(std::ostream &out) const
void SetTime(double t)
Set physical time (for time-dependent simulations)
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
The classical explicit forth-order Runge-Kutta method, RK4.
void velocity_function(const Vector &x, Vector &v)
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Third-order, strong stability preserving (SSP) Runge-Kutta method.
NURBSExtension * NURBSext
void PrintOptions(std::ostream &out) const
void ProjectCoefficient(Coefficient &coeff)
int open(const char hostname[], int port)
double u0_function(const Vector &x)
class for C-function coefficient
The classical forward Euler method.
double inflow_function(const Vector &x)
Arbitrary order "L2-conforming" discontinuous finite elements.