61 int main(
int argc,
char *argv[])
64 MPI_Session mpi(argc, argv);
68 const char *mesh_file =
"../data/periodic-square.mesh";
69 int ser_ref_levels = 0;
70 int par_ref_levels = 1;
72 int ode_solver_type = 4;
76 bool visualization =
true;
80 cout.precision(precision);
82 OptionsParser args(argc, argv);
83 args.AddOption(&mesh_file,
"-m",
"--mesh",
85 args.AddOption(&
problem,
"-p",
"--problem",
86 "Problem setup to use. See options in velocity_function().");
87 args.AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
88 "Number of times to refine the mesh uniformly before parallel"
89 " partitioning, -1 for auto.");
90 args.AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
91 "Number of times to refine the mesh uniformly after parallel"
93 args.AddOption(&order,
"-o",
"--order",
94 "Order (degree) of the finite elements.");
95 args.AddOption(&ode_solver_type,
"-s",
"--ode-solver",
96 "ODE solver: 1 - Forward Euler,\n\t"
97 " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6.");
98 args.AddOption(&t_final,
"-tf",
"--t-final",
99 "Final time; start time is 0.");
100 args.AddOption(&dt,
"-dt",
"--time-step",
101 "Time step. Positive number skips CFL timestep calculation.");
102 args.AddOption(&cfl,
"-c",
"--cfl-number",
103 "CFL number for timestep calculation.");
104 args.AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
105 "--no-visualization",
106 "Enable or disable GLVis visualization.");
107 args.AddOption(&vis_steps,
"-vs",
"--visualization-steps",
108 "Visualize every n-th timestep.");
113 if (mpi.Root()) { args.PrintUsage(cout); }
116 if (mpi.Root()) { args.PrintOptions(cout); }
120 Mesh mesh(mesh_file, 1, 1);
121 const int dim = mesh.Dimension();
123 MFEM_ASSERT(dim == 2,
"Need a two-dimensional mesh for the problem definition");
127 ODESolver *ode_solver = NULL;
128 switch (ode_solver_type)
130 case 1: ode_solver =
new ForwardEulerSolver;
break;
131 case 2: ode_solver =
new RK2Solver(1.0);
break;
132 case 3: ode_solver =
new RK3SSPSolver;
break;
133 case 4: ode_solver =
new RK4Solver;
break;
134 case 6: ode_solver =
new RK6Solver;
break;
138 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
146 for (
int lev = 0; lev < ser_ref_levels; lev++)
148 mesh.UniformRefinement();
154 ParMesh pmesh(MPI_COMM_WORLD, mesh);
156 for (
int lev = 0; lev < par_ref_levels; lev++)
158 pmesh.UniformRefinement();
165 ParFiniteElementSpace fes(&pmesh, &fec);
167 ParFiniteElementSpace dfes(&pmesh, &fec, dim, Ordering::byNODES);
169 ParFiniteElementSpace vfes(&pmesh, &fec,
num_equation, Ordering::byNODES);
172 MFEM_ASSERT(fes.GetOrdering() == Ordering::byNODES,
"");
175 if (mpi.Root()) { cout <<
"Number of unknowns: " << glob_size << endl; }
183 for (
int k = 0; k <=
num_equation; k++) { offsets[k] = k * vfes.GetNDofs(); }
184 BlockVector u_block(offsets);
187 ParGridFunction mom(&dfes, u_block.GetData() + offsets[1]);
191 ParGridFunction sol(&vfes, u_block.GetData());
192 sol.ProjectCoefficient(u0);
196 ostringstream mesh_name;
197 mesh_name <<
"vortex-mesh." << setfill(
'0') << setw(6) << mpi.WorldRank();
198 ofstream mesh_ofs(mesh_name.str().c_str());
199 mesh_ofs.precision(precision);
204 ParGridFunction uk(&fes, u_block.GetBlock(k));
205 ostringstream sol_name;
206 sol_name <<
"vortex-" << k <<
"-init."
207 << setfill(
'0') << setw(6) << mpi.WorldRank();
208 ofstream sol_ofs(sol_name.str().c_str());
209 sol_ofs.precision(precision);
216 MixedBilinearForm Aflux(&dfes, &fes);
220 ParNonlinearForm A(&vfes);
236 MPI_Barrier(pmesh.GetComm());
237 sout.open(vishost, visport);
242 cout <<
"Unable to connect to GLVis server at "
243 << vishost <<
':' << visport << endl;
245 visualization =
false;
246 if (mpi.Root()) { cout <<
"GLVis visualization disabled.\n"; }
250 sout <<
"parallel " << mpi.WorldSize() <<
" " << mpi.WorldRank() <<
"\n";
251 sout.precision(precision);
252 sout <<
"solution\n" << pmesh << mom;
257 cout <<
"GLVis visualization paused."
258 <<
" Press space (in the GLVis window) to resume it.\n";
267 double my_hmin = pmesh.GetElementSize(0, 1);
268 for (
int i = 1; i < pmesh.GetNE(); i++)
270 my_hmin = min(pmesh.GetElementSize(i, 1), my_hmin);
273 MPI_Allreduce(&my_hmin, &hmin, 1, MPI_DOUBLE, MPI_MIN, pmesh.GetComm());
282 ode_solver->Init(euler);
289 Vector z(sol.Size());
293 double all_max_char_speed;
295 1, MPI_DOUBLE, MPI_MAX, pmesh.GetComm());
303 for (
int ti = 0; !done; )
305 double dt_real = min(dt, t_final - t);
307 ode_solver->Step(sol, t, dt_real);
312 double all_max_char_speed;
314 1, MPI_DOUBLE, MPI_MAX, pmesh.GetComm());
321 done = (t >= t_final - 1e-8*dt);
322 if (done || ti % vis_steps == 0)
326 cout <<
"time step: " << ti <<
", time: " << t << endl;
330 MPI_Barrier(pmesh.GetComm());
331 sout <<
"parallel " << mpi.WorldSize() <<
" " << mpi.WorldRank() <<
"\n";
332 sout <<
"solution\n" << pmesh << mom << flush;
338 if (mpi.Root()) { cout <<
" done, " <<
tic_toc.
RealTime() <<
"s." << endl; }
344 ParGridFunction uk(&fes, u_block.GetBlock(k));
345 ostringstream sol_name;
346 sol_name <<
"vortex-" << k <<
"-final."
347 << setfill(
'0') << setw(6) << mpi.WorldRank();
348 ofstream sol_ofs(sol_name.str().c_str());
349 sol_ofs.precision(precision);
356 const double error = sol.ComputeLpError(2, u0);
357 if (mpi.Root()) { cout <<
"Solution error: " << error << endl; }
void Stop()
Stop the stopwatch.
const double specific_heat_ratio
void Start()
Clear the elapsed time and start the stopwatch.
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
void InitialCondition(const Vector &x, Vector &y)
const double gas_constant
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.