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;
95 cout.precision(precision);
98 args.
AddOption(&mesh_file,
"-m",
"--mesh",
101 "Problem setup to use. See options in velocity_function().");
102 args.
AddOption(&ref_levels,
"-r",
"--refine",
103 "Number of times to refine the mesh uniformly.");
105 "Order (degree) of the finite elements.");
106 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
107 "ODE solver: 1 - Forward Euler,\n\t"
108 " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6.");
109 args.
AddOption(&t_final,
"-tf",
"--t-final",
110 "Final time; start time is 0.");
111 args.
AddOption(&dt,
"-dt",
"--time-step",
113 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
114 "--no-visualization",
115 "Enable or disable GLVis visualization.");
116 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
117 "--no-visit-datafiles",
118 "Save data files for VisIt (visit.llnl.gov) visualization.");
119 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
121 "Use binary (Sidre) or ascii format for VisIt data files.");
122 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
123 "Visualize every n-th timestep.");
134 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
140 switch (ode_solver_type)
143 case 2: ode_solver =
new RK2Solver(1.0);
break;
145 case 4: ode_solver =
new RK4Solver;
break;
146 case 6: ode_solver =
new RK6Solver;
break;
148 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
156 for (
int lev = 0; lev < ref_levels; lev++)
171 cout <<
"Number of unknowns: " << fes.
GetVSize() << endl;
207 ofstream omesh(
"ex9.mesh");
208 omesh.precision(precision);
210 ofstream osol(
"ex9-init.gf");
211 osol.precision(precision);
222 #ifdef MFEM_USE_SIDRE
225 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
242 char vishost[] =
"localhost";
244 sout.
open(vishost, visport);
247 cout <<
"Unable to connect to GLVis server at "
248 << vishost <<
':' << visport << endl;
249 visualization =
false;
250 cout <<
"GLVis visualization disabled.\n";
254 sout.precision(precision);
255 sout <<
"solution\n" << *mesh << u;
258 cout <<
"GLVis visualization paused."
259 <<
" Press space (in the GLVis window) to resume it.\n";
270 ode_solver->
Init(adv);
273 for (
int ti = 0; !done; )
275 double dt_real = min(dt, t_final - t);
276 ode_solver->
Step(u, t, dt_real);
279 done = (t >= t_final - 1e-8*dt);
281 if (done || ti % vis_steps == 0)
283 cout <<
"time step: " << ti <<
", time: " << t << endl;
287 sout <<
"solution\n" << *mesh << u << flush;
302 ofstream osol(
"ex9-final.gf");
303 osol.precision(precision);
319 M_solver.SetPreconditioner(M_prec);
320 M_solver.SetOperator(M);
322 M_solver.iterative_mode =
false;
323 M_solver.SetRelTol(1e-9);
324 M_solver.SetAbsTol(0.0);
325 M_solver.SetMaxIter(100);
326 M_solver.SetPrintLevel(0);
345 for (
int i = 0; i <
dim; i++)
358 case 1: v(0) = 1.0;
break;
359 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
360 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
369 const double w = M_PI/2;
372 case 1: v(0) = 1.0;
break;
373 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
374 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
381 const double w = M_PI/2;
382 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
386 case 1: v(0) = 1.0;
break;
387 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
388 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
402 for (
int i = 0; i <
dim; i++)
416 return exp(-40.*pow(X(0)-0.5,2));
420 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
423 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
427 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
428 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
434 double x_ = X(0), y_ = X(1), rho, phi;
437 return pow(sin(M_PI*rho),2)*sin(3*phi);
441 const double f = M_PI;
442 return sin(f*X(0))*sin(f*X(1));
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
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 Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Base abstract class for time dependent operators.
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
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
int Size() const
Returns the size of the vector.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
int main(int argc, char *argv[])
Data collection with Sidre routines following the Conduit mesh blueprint specification.
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)
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
Optional NURBS mesh extension.
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
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
The classical forward Euler method.
double inflow_function(const Vector &x)
Arbitrary order "L2-conforming" discontinuous finite elements.
virtual void Print(std::ostream &out=std::cout) const