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