78 int main(
int argc,
char *argv[])
82 MPI_Init(&argc, &argv);
83 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
84 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
88 const char *mesh_file =
"../data/periodic-hexagon.mesh";
89 int ser_ref_levels = 2;
90 int par_ref_levels = 0;
92 int ode_solver_type = 4;
93 double t_final = 10.0;
95 bool visualization =
true;
101 cout.precision(precision);
104 args.
AddOption(&mesh_file,
"-m",
"--mesh",
105 "Mesh file to use.");
107 "Problem setup to use. See options in velocity_function().");
108 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
109 "Number of times to refine the mesh uniformly in serial.");
110 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
111 "Number of times to refine the mesh uniformly in parallel.");
113 "Order (degree) of the finite elements.");
114 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
115 "ODE solver: 1 - Forward Euler,\n\t"
116 " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6.");
117 args.
AddOption(&t_final,
"-tf",
"--t-final",
118 "Final time; start time is 0.");
119 args.
AddOption(&dt,
"-dt",
"--time-step",
121 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
122 "--no-visualization",
123 "Enable or disable GLVis visualization.");
124 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
125 "--no-visit-datafiles",
126 "Save data files for VisIt (visit.llnl.gov) visualization.");
127 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
129 "Use binary (Sidre) or ascii format for VisIt data files.");
130 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
131 "Visualize every n-th timestep.");
149 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
155 switch (ode_solver_type)
158 case 2: ode_solver =
new RK2Solver(1.0);
break;
160 case 4: ode_solver =
new RK4Solver;
break;
161 case 6: ode_solver =
new RK6Solver;
break;
165 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
176 for (
int lev = 0; lev < ser_ref_levels; lev++)
191 for (
int lev = 0; lev < par_ref_levels; lev++)
204 cout <<
"Number of unknowns: " << global_vSize << endl;
246 ostringstream mesh_name, sol_name;
247 mesh_name <<
"ex9-mesh." << setfill(
'0') << setw(6) << myid;
248 sol_name <<
"ex9-init." << setfill(
'0') << setw(6) << myid;
249 ofstream omesh(mesh_name.str().c_str());
250 omesh.precision(precision);
252 ofstream osol(sol_name.str().c_str());
253 osol.precision(precision);
264 #ifdef MFEM_USE_SIDRE
267 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
286 char vishost[] =
"localhost";
288 sout.
open(vishost, visport);
292 cout <<
"Unable to connect to GLVis server at "
293 << vishost <<
':' << visport << endl;
294 visualization =
false;
297 cout <<
"GLVis visualization disabled.\n";
302 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
303 sout.precision(precision);
304 sout <<
"solution\n" << *pmesh << *u;
308 cout <<
"GLVis visualization paused."
309 <<
" Press space (in the GLVis window) to resume it.\n";
320 ode_solver->
Init(adv);
323 for (
int ti = 0; !done; )
325 double dt_real = min(dt, t_final - t);
326 ode_solver->
Step(*U, t, dt_real);
329 done = (t >= t_final - 1e-8*dt);
331 if (done || ti % vis_steps == 0)
335 cout <<
"time step: " << ti <<
", time: " << t << endl;
344 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
345 sout <<
"solution\n" << *pmesh << *u << flush;
361 ostringstream sol_name;
362 sol_name <<
"ex9-final." << setfill(
'0') << setw(6) << myid;
363 ofstream osol(sol_name.str().c_str());
364 osol.precision(precision);
391 M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height())
393 M_prec.SetType(HypreSmoother::Jacobi);
420 for (
int i = 0; i <
dim; i++)
433 case 1: v(0) = 1.0;
break;
434 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
435 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
444 const double w = M_PI/2;
447 case 1: v(0) = 1.0;
break;
448 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
449 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
456 const double w = M_PI/2;
457 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
461 case 1: v(0) = 1.0;
break;
462 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
463 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
477 for (
int i = 0; i <
dim; i++)
491 return exp(-40.*pow(X(0)-0.5,2));
495 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
498 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
502 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
503 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
509 double x_ = X(0), y_ = X(1), rho, phi;
512 return pow(sin(M_PI*rho),2)*sin(3*phi);
516 const double f = M_PI;
517 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.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
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]...
virtual void SetTime(const double _t)
Set the current time.
int Size() const
Returns the size of the vector.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
virtual void Save(std::ostream &out) const
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
int main(int argc, char *argv[])
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Data collection with Sidre routines following the Conduit mesh blueprint specification.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetPrintLevel(int print_lvl)
virtual void Save()
Save the collection to disk.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
HYPRE_Int GlobalTrueVSize() const
virtual void Print(std::ostream &out=mfem::out) const
Parallel smoothers in hypre.
void PrintUsage(std::ostream &out) const
void SetTime(double t)
Set physical time (for time-dependent simulations)
Wrapper for hypre's parallel vector class.
The classical explicit forth-order Runge-Kutta method, RK4.
void SetAbsTol(double atol)
void SetRelTol(double rtol)
void velocity_function(const Vector &x, Vector &v)
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
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.
virtual void Mult(const Vector &x, Vector &y) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void PrintOptions(std::ostream &out) const
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
int open(const char hostname[], int port)
double u0_function(const Vector &x)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
class for C-function coefficient
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Class for parallel grid function.
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
double inflow_function(const Vector &x)
Arbitrary order "L2-conforming" discontinuous finite elements.