MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex16p.cpp
Go to the documentation of this file.
1 // MFEM Example 16 - Parallel Version
2 //
3 // Compile with: make ex16p
4 //
5 // Sample runs: mpirun -np 4 ex16p
6 // mpirun -np 4 ex16p -m ../data/inline-tri.mesh
7 // mpirun -np 4 ex16p -m ../data/disc-nurbs.mesh -tf 2
8 // mpirun -np 4 ex16p -s 1 -a 0.0 -k 1.0
9 // mpirun -np 4 ex16p -s 2 -a 1.0 -k 0.0
10 // mpirun -np 8 ex16p -s 3 -a 0.5 -k 0.5 -o 4
11 // mpirun -np 4 ex16p -s 14 -dt 1.0e-4 -tf 4.0e-2 -vs 40
12 // mpirun -np 16 ex16p -m ../data/fichera-q2.mesh
13 // mpirun -np 16 ex16p -m ../data/fichera-mixed.mesh
14 // mpirun -np 16 ex16p -m ../data/escher-p2.mesh
15 // mpirun -np 8 ex16p -m ../data/beam-tet.mesh -tf 10 -dt 0.1
16 // mpirun -np 4 ex16p -m ../data/amr-quad.mesh -o 4 -rs 0 -rp 0
17 // mpirun -np 4 ex16p -m ../data/amr-hex.mesh -o 2 -rs 0 -rp 0
18 //
19 // Description: This example solves a time dependent nonlinear heat equation
20 // problem of the form du/dt = C(u), with a non-linear diffusion
21 // operator C(u) = \nabla \cdot (\kappa + \alpha u) \nabla u.
22 //
23 // The example demonstrates the use of nonlinear operators (the
24 // class ConductionOperator defining C(u)), as well as their
25 // implicit time integration. Note that implementing the method
26 // ConductionOperator::ImplicitSolve is the only requirement for
27 // high-order implicit (SDIRK) time integration. Optional saving
28 // with ADIOS2 (adios2.readthedocs.io) is also illustrated.
29 //
30 // We recommend viewing examples 2, 9 and 10 before viewing this
31 // example.
32 
33 #include "mfem.hpp"
34 #include <fstream>
35 #include <iostream>
36 
37 using namespace std;
38 using namespace mfem;
39 
40 /** After spatial discretization, the conduction model can be written as:
41  *
42  * du/dt = M^{-1}(-Ku)
43  *
44  * where u is the vector representing the temperature, M is the mass matrix,
45  * and K is the diffusion operator with diffusivity depending on u:
46  * (\kappa + \alpha u).
47  *
48  * Class ConductionOperator represents the right-hand side of the above ODE.
49  */
50 class ConductionOperator : public TimeDependentOperator
51 {
52 protected:
53  ParFiniteElementSpace &fespace;
54  Array<int> ess_tdof_list; // this list remains empty for pure Neumann b.c.
55 
56  ParBilinearForm *M;
57  ParBilinearForm *K;
58 
59  HypreParMatrix Mmat;
60  HypreParMatrix Kmat;
61  HypreParMatrix *T; // T = M + dt K
62  double current_dt;
63 
64  CGSolver M_solver; // Krylov solver for inverting the mass matrix M
65  HypreSmoother M_prec; // Preconditioner for the mass matrix M
66 
67  CGSolver T_solver; // Implicit solver for T = M + dt K
68  HypreSmoother T_prec; // Preconditioner for the implicit solver
69 
70  double alpha, kappa;
71 
72  mutable Vector z; // auxiliary vector
73 
74 public:
75  ConductionOperator(ParFiniteElementSpace &f, double alpha, double kappa,
76  const Vector &u);
77 
78  virtual void Mult(const Vector &u, Vector &du_dt) const;
79  /** Solve the Backward-Euler equation: k = f(u + dt*k, t), for the unknown k.
80  This is the only requirement for high-order SDIRK implicit integration.*/
81  virtual void ImplicitSolve(const double dt, const Vector &u, Vector &k);
82 
83  /// Update the diffusion BilinearForm K using the given true-dof vector `u`.
84  void SetParameters(const Vector &u);
85 
86  virtual ~ConductionOperator();
87 };
88 
89 double InitialTemperature(const Vector &x);
90 
91 int main(int argc, char *argv[])
92 {
93  // 1. Initialize MPI.
94  int num_procs, myid;
95  MPI_Init(&argc, &argv);
96  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
97  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
98 
99  // 2. Parse command-line options.
100  const char *mesh_file = "../data/star.mesh";
101  int ser_ref_levels = 2;
102  int par_ref_levels = 1;
103  int order = 2;
104  int ode_solver_type = 3;
105  double t_final = 0.5;
106  double dt = 1.0e-2;
107  double alpha = 1.0e-2;
108  double kappa = 0.5;
109  bool visualization = true;
110  bool visit = false;
111  int vis_steps = 5;
112  bool adios2 = false;
113 
114  int precision = 8;
115  cout.precision(precision);
116 
117  OptionsParser args(argc, argv);
118  args.AddOption(&mesh_file, "-m", "--mesh",
119  "Mesh file to use.");
120  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
121  "Number of times to refine the mesh uniformly in serial.");
122  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
123  "Number of times to refine the mesh uniformly in parallel.");
124  args.AddOption(&order, "-o", "--order",
125  "Order (degree) of the finite elements.");
126  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
127  "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
128  "\t 11 - Forward Euler, 12 - RK2, 13 - RK3 SSP, 14 - RK4.");
129  args.AddOption(&t_final, "-tf", "--t-final",
130  "Final time; start time is 0.");
131  args.AddOption(&dt, "-dt", "--time-step",
132  "Time step.");
133  args.AddOption(&alpha, "-a", "--alpha",
134  "Alpha coefficient.");
135  args.AddOption(&kappa, "-k", "--kappa",
136  "Kappa coefficient offset.");
137  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
138  "--no-visualization",
139  "Enable or disable GLVis visualization.");
140  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
141  "--no-visit-datafiles",
142  "Save data files for VisIt (visit.llnl.gov) visualization.");
143  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
144  "Visualize every n-th timestep.");
145  args.AddOption(&adios2, "-adios2", "--adios2-streams", "-no-adios2",
146  "--no-adios2-streams",
147  "Save data using adios2 streams.");
148  args.Parse();
149  if (!args.Good())
150  {
151  args.PrintUsage(cout);
152  MPI_Finalize();
153  return 1;
154  }
155 
156  if (myid == 0)
157  {
158  args.PrintOptions(cout);
159  }
160 
161  // 3. Read the serial mesh from the given mesh file on all processors. We can
162  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
163  // with the same code.
164  Mesh *mesh = new Mesh(mesh_file, 1, 1);
165  int dim = mesh->Dimension();
166 
167  // 4. Define the ODE solver used for time integration. Several implicit
168  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
169  // explicit Runge-Kutta methods are available.
170  ODESolver *ode_solver;
171  switch (ode_solver_type)
172  {
173  // Implicit L-stable methods
174  case 1: ode_solver = new BackwardEulerSolver; break;
175  case 2: ode_solver = new SDIRK23Solver(2); break;
176  case 3: ode_solver = new SDIRK33Solver; break;
177  // Explicit methods
178  case 11: ode_solver = new ForwardEulerSolver; break;
179  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
180  case 13: ode_solver = new RK3SSPSolver; break;
181  case 14: ode_solver = new RK4Solver; break;
182  case 15: ode_solver = new GeneralizedAlphaSolver(0.5); break;
183  // Implicit A-stable methods (not L-stable)
184  case 22: ode_solver = new ImplicitMidpointSolver; break;
185  case 23: ode_solver = new SDIRK23Solver; break;
186  case 24: ode_solver = new SDIRK34Solver; break;
187  default:
188  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
189  delete mesh;
190  return 3;
191  }
192 
193  // 5. Refine the mesh in serial to increase the resolution. In this example
194  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
195  // a command-line parameter.
196  for (int lev = 0; lev < ser_ref_levels; lev++)
197  {
198  mesh->UniformRefinement();
199  }
200 
201  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
202  // this mesh further in parallel to increase the resolution. Once the
203  // parallel mesh is defined, the serial mesh can be deleted.
204  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
205  delete mesh;
206  for (int lev = 0; lev < par_ref_levels; lev++)
207  {
208  pmesh->UniformRefinement();
209  }
210 
211  // 7. Define the vector finite element space representing the current and the
212  // initial temperature, u_ref.
213  H1_FECollection fe_coll(order, dim);
214  ParFiniteElementSpace fespace(pmesh, &fe_coll);
215 
216  HYPRE_BigInt fe_size = fespace.GlobalTrueVSize();
217  if (myid == 0)
218  {
219  cout << "Number of temperature unknowns: " << fe_size << endl;
220  }
221 
222  ParGridFunction u_gf(&fespace);
223 
224  // 8. Set the initial conditions for u. All boundaries are considered
225  // natural.
227  u_gf.ProjectCoefficient(u_0);
228  Vector u;
229  u_gf.GetTrueDofs(u);
230 
231  // 9. Initialize the conduction operator and the VisIt visualization.
232  ConductionOperator oper(fespace, alpha, kappa, u);
233 
234  u_gf.SetFromTrueDofs(u);
235  {
236  ostringstream mesh_name, sol_name;
237  mesh_name << "ex16-mesh." << setfill('0') << setw(6) << myid;
238  sol_name << "ex16-init." << setfill('0') << setw(6) << myid;
239  ofstream omesh(mesh_name.str().c_str());
240  omesh.precision(precision);
241  pmesh->Print(omesh);
242  ofstream osol(sol_name.str().c_str());
243  osol.precision(precision);
244  u_gf.Save(osol);
245  }
246 
247  VisItDataCollection visit_dc("Example16-Parallel", pmesh);
248  visit_dc.RegisterField("temperature", &u_gf);
249  if (visit)
250  {
251  visit_dc.SetCycle(0);
252  visit_dc.SetTime(0.0);
253  visit_dc.Save();
254  }
255 
256  // Optionally output a BP (binary pack) file using ADIOS2. This can be
257  // visualized with the ParaView VTX reader.
258 #ifdef MFEM_USE_ADIOS2
259  ADIOS2DataCollection* adios2_dc = NULL;
260  if (adios2)
261  {
262  std::string postfix(mesh_file);
263  postfix.erase(0, std::string("../data/").size() );
264  postfix += "_o" + std::to_string(order);
265  postfix += "_solver" + std::to_string(ode_solver_type);
266  const std::string collection_name = "ex16-p-" + postfix + ".bp";
267 
268  adios2_dc = new ADIOS2DataCollection(MPI_COMM_WORLD, collection_name, pmesh);
269  adios2_dc->SetParameter("SubStreams", std::to_string(num_procs/2) );
270  adios2_dc->RegisterField("temperature", &u_gf);
271  adios2_dc->SetCycle(0);
272  adios2_dc->SetTime(0.0);
273  adios2_dc->Save();
274  }
275 #endif
276 
277  socketstream sout;
278  if (visualization)
279  {
280  char vishost[] = "localhost";
281  int visport = 19916;
282  sout.open(vishost, visport);
283  sout << "parallel " << num_procs << " " << myid << endl;
284  int good = sout.good(), all_good;
285  MPI_Allreduce(&good, &all_good, 1, MPI_INT, MPI_MIN, pmesh->GetComm());
286  if (!all_good)
287  {
288  sout.close();
289  visualization = false;
290  if (myid == 0)
291  {
292  cout << "Unable to connect to GLVis server at "
293  << vishost << ':' << visport << endl;
294  cout << "GLVis visualization disabled.\n";
295  }
296  }
297  else
298  {
299  sout.precision(precision);
300  sout << "solution\n" << *pmesh << u_gf;
301  sout << "pause\n";
302  sout << flush;
303  if (myid == 0)
304  {
305  cout << "GLVis visualization paused."
306  << " Press space (in the GLVis window) to resume it.\n";
307  }
308  }
309  }
310 
311  // 10. Perform time-integration (looping over the time iterations, ti, with a
312  // time-step dt).
313  ode_solver->Init(oper);
314  double t = 0.0;
315 
316  bool last_step = false;
317  for (int ti = 1; !last_step; ti++)
318  {
319  if (t + dt >= t_final - dt/2)
320  {
321  last_step = true;
322  }
323 
324  ode_solver->Step(u, t, dt);
325 
326  if (last_step || (ti % vis_steps) == 0)
327  {
328  if (myid == 0)
329  {
330  cout << "step " << ti << ", t = " << t << endl;
331  }
332 
333  u_gf.SetFromTrueDofs(u);
334  if (visualization)
335  {
336  sout << "parallel " << num_procs << " " << myid << "\n";
337  sout << "solution\n" << *pmesh << u_gf << flush;
338  }
339 
340  if (visit)
341  {
342  visit_dc.SetCycle(ti);
343  visit_dc.SetTime(t);
344  visit_dc.Save();
345  }
346 
347 #ifdef MFEM_USE_ADIOS2
348  if (adios2)
349  {
350  adios2_dc->SetCycle(ti);
351  adios2_dc->SetTime(t);
352  adios2_dc->Save();
353  }
354 #endif
355  }
356  oper.SetParameters(u);
357  }
358 
359 #ifdef MFEM_USE_ADIOS2
360  if (adios2)
361  {
362  delete adios2_dc;
363  }
364 #endif
365 
366  // 11. Save the final solution in parallel. This output can be viewed later
367  // using GLVis: "glvis -np <np> -m ex16-mesh -g ex16-final".
368  {
369  ostringstream sol_name;
370  sol_name << "ex16-final." << setfill('0') << setw(6) << myid;
371  ofstream osol(sol_name.str().c_str());
372  osol.precision(precision);
373  u_gf.Save(osol);
374  }
375 
376  // 12. Free the used memory.
377  delete ode_solver;
378  delete pmesh;
379 
380  MPI_Finalize();
381 
382  return 0;
383 }
384 
385 ConductionOperator::ConductionOperator(ParFiniteElementSpace &f, double al,
386  double kap, const Vector &u)
387  : TimeDependentOperator(f.GetTrueVSize(), 0.0), fespace(f), M(NULL), K(NULL),
388  T(NULL), current_dt(0.0),
389  M_solver(f.GetComm()), T_solver(f.GetComm()), z(height)
390 {
391  const double rel_tol = 1e-8;
392 
393  M = new ParBilinearForm(&fespace);
394  M->AddDomainIntegrator(new MassIntegrator());
395  M->Assemble(0); // keep sparsity pattern of M and K the same
396  M->FormSystemMatrix(ess_tdof_list, Mmat);
397 
398  M_solver.iterative_mode = false;
399  M_solver.SetRelTol(rel_tol);
400  M_solver.SetAbsTol(0.0);
401  M_solver.SetMaxIter(100);
402  M_solver.SetPrintLevel(0);
403  M_prec.SetType(HypreSmoother::Jacobi);
404  M_solver.SetPreconditioner(M_prec);
405  M_solver.SetOperator(Mmat);
406 
407  alpha = al;
408  kappa = kap;
409 
410  T_solver.iterative_mode = false;
411  T_solver.SetRelTol(rel_tol);
412  T_solver.SetAbsTol(0.0);
413  T_solver.SetMaxIter(100);
414  T_solver.SetPrintLevel(0);
415  T_solver.SetPreconditioner(T_prec);
416 
417  SetParameters(u);
418 }
419 
420 void ConductionOperator::Mult(const Vector &u, Vector &du_dt) const
421 {
422  // Compute:
423  // du_dt = M^{-1}*-K(u)
424  // for du_dt
425  Kmat.Mult(u, z);
426  z.Neg(); // z = -z
427  M_solver.Mult(z, du_dt);
428 }
429 
430 void ConductionOperator::ImplicitSolve(const double dt,
431  const Vector &u, Vector &du_dt)
432 {
433  // Solve the equation:
434  // du_dt = M^{-1}*[-K(u + dt*du_dt)]
435  // for du_dt
436  if (!T)
437  {
438  T = Add(1.0, Mmat, dt, Kmat);
439  current_dt = dt;
440  T_solver.SetOperator(*T);
441  }
442  MFEM_VERIFY(dt == current_dt, ""); // SDIRK methods use the same dt
443  Kmat.Mult(u, z);
444  z.Neg();
445  T_solver.Mult(z, du_dt);
446 }
447 
448 void ConductionOperator::SetParameters(const Vector &u)
449 {
450  ParGridFunction u_alpha_gf(&fespace);
451  u_alpha_gf.SetFromTrueDofs(u);
452  for (int i = 0; i < u_alpha_gf.Size(); i++)
453  {
454  u_alpha_gf(i) = kappa + alpha*u_alpha_gf(i);
455  }
456 
457  delete K;
458  K = new ParBilinearForm(&fespace);
459 
460  GridFunctionCoefficient u_coeff(&u_alpha_gf);
461 
462  K->AddDomainIntegrator(new DiffusionIntegrator(u_coeff));
463  K->Assemble(0); // keep sparsity pattern of M and K the same
464  K->FormSystemMatrix(ess_tdof_list, Kmat);
465  delete T;
466  T = NULL; // re-compute T on the next ImplicitSolve
467 }
468 
469 ConductionOperator::~ConductionOperator()
470 {
471  delete T;
472  delete M;
473  delete K;
474 }
475 
476 double InitialTemperature(const Vector &x)
477 {
478  if (x.Norml2() < 0.5)
479  {
480  return 2.0;
481  }
482  else
483  {
484  return 1.0;
485  }
486 }
Conjugate gradient method.
Definition: solvers.hpp:316
double InitialTemperature(const Vector &x)
Definition: ex16.cpp:382
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Base abstract class for first order time dependent operators.
Definition: operator.hpp:282
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:783
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]...
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:493
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:389
int close()
Close the socketstream.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1931
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:164
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
constexpr char vishost[]
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
Data collection with VisIt I/O routines.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
void SetParameter(const std::string key, const std::string value) noexcept
Parallel smoothers in hypre.
Definition: hypre.hpp:840
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void SetTime(double t)
Set physical time (for time-dependent simulations)
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
int precision
Precision (number of digits) used for the text output of doubles.
MPI_Comm GetComm() const
Definition: pmesh.hpp:276
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)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:149
HYPRE_Int HYPRE_BigInt
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:402
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
Class for parallel bilinear form.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
RefCoord t[3]
const double alpha
Definition: ex15.cpp:369
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The classical forward Euler method.
Definition: ode.hpp:116
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
int main()
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150