MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex9.cpp
Go to the documentation of this file.
1 // MFEM Example 9
2 //
3 // Compile with: make ex9
4 //
5 // Sample runs:
6 // ex9 -m ../data/periodic-segment.mesh -p 0 -r 2 -dt 0.005
7 // ex9 -m ../data/periodic-square.mesh -p 0 -r 2 -dt 0.01 -tf 10
8 // ex9 -m ../data/periodic-hexagon.mesh -p 0 -r 2 -dt 0.01 -tf 10
9 // ex9 -m ../data/periodic-square.mesh -p 1 -r 2 -dt 0.005 -tf 9
10 // ex9 -m ../data/periodic-hexagon.mesh -p 1 -r 2 -dt 0.005 -tf 9
11 // ex9 -m ../data/amr-quad.mesh -p 1 -r 2 -dt 0.002 -tf 9
12 // ex9 -m ../data/star-q3.mesh -p 1 -r 2 -dt 0.005 -tf 9
13 // ex9 -m ../data/disc-nurbs.mesh -p 1 -r 3 -dt 0.005 -tf 9
14 // ex9 -m ../data/disc-nurbs.mesh -p 2 -r 3 -dt 0.005 -tf 9
15 // ex9 -m ../data/periodic-square.mesh -p 3 -r 4 -dt 0.0025 -tf 9 -vs 20
16 // ex9 -m ../data/periodic-cube.mesh -p 0 -r 2 -o 2 -dt 0.02 -tf 8
17 //
18 // Description: This example code solves the time-dependent advection equation
19 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
20 // u0(x)=u(0,x) is a given initial condition.
21 //
22 // The example demonstrates the use of Discontinuous Galerkin (DG)
23 // bilinear forms in MFEM (face integrators), the use of explicit
24 // ODE time integrators, the definition of periodic boundary
25 // conditions through periodic meshes, as well as the use of GLVis
26 // for persistent visualization of a time-evolving solution. The
27 // saving of time-dependent data files for external visualization
28 // with VisIt (visit.llnl.gov) is also illustrated.
29 
30 #include "mfem.hpp"
31 #include <fstream>
32 #include <iostream>
33 #include <algorithm>
34 
35 using namespace std;
36 using namespace mfem;
37 
38 // Choice for the problem setup. The fluid velocity, initial condition and
39 // inflow boundary condition are chosen based on this parameter.
40 int problem;
41 
42 // Velocity coefficient
43 void velocity_function(const Vector &x, Vector &v);
44 
45 // Initial condition
46 double u0_function(const Vector &x);
47 
48 // Inflow boundary condition
49 double inflow_function(const Vector &x);
50 
51 // Mesh bounding box
53 
54 
55 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
56  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
57  and advection matrices, and b describes the flow on the boundary. This can
58  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
59  used to evaluate the right-hand side. */
60 class FE_Evolution : public TimeDependentOperator
61 {
62 private:
63  SparseMatrix &M, &K;
64  const Vector &b;
65  DSmoother M_prec;
66  CGSolver M_solver;
67 
68  mutable Vector z;
69 
70 public:
71  FE_Evolution(SparseMatrix &_M, SparseMatrix &_K, const Vector &_b);
72 
73  virtual void Mult(const Vector &x, Vector &y) const;
74 
75  virtual ~FE_Evolution() { }
76 };
77 
78 
79 int main(int argc, char *argv[])
80 {
81  // 1. Parse command-line options.
82  problem = 0;
83  const char *mesh_file = "../data/periodic-hexagon.mesh";
84  int ref_levels = 2;
85  int order = 3;
86  int ode_solver_type = 4;
87  double t_final = 10.0;
88  double dt = 0.01;
89  bool visualization = true;
90  bool visit = false;
91  bool binary = false;
92  int vis_steps = 5;
93 
94  int precision = 8;
95  cout.precision(precision);
96 
97  OptionsParser args(argc, argv);
98  args.AddOption(&mesh_file, "-m", "--mesh",
99  "Mesh file to use.");
100  args.AddOption(&problem, "-p", "--problem",
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.");
104  args.AddOption(&order, "-o", "--order",
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",
112  "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",
120  "--ascii-datafiles",
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.");
124  args.Parse();
125  if (!args.Good())
126  {
127  args.PrintUsage(cout);
128  return 1;
129  }
130  args.PrintOptions(cout);
131 
132  // 2. Read the mesh from the given mesh file. We can handle geometrically
133  // periodic meshes in this code.
134  Mesh *mesh = new Mesh(mesh_file, 1, 1);
135  int dim = mesh->Dimension();
136 
137  // 3. Define the ODE solver used for time integration. Several explicit
138  // Runge-Kutta methods are available.
139  ODESolver *ode_solver = NULL;
140  switch (ode_solver_type)
141  {
142  case 1: ode_solver = new ForwardEulerSolver; break;
143  case 2: ode_solver = new RK2Solver(1.0); break;
144  case 3: ode_solver = new RK3SSPSolver; break;
145  case 4: ode_solver = new RK4Solver; break;
146  case 6: ode_solver = new RK6Solver; break;
147  default:
148  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
149  return 3;
150  }
151 
152  // 4. Refine the mesh to increase the resolution. In this example we do
153  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
154  // command-line parameter. If the mesh is of NURBS type, we convert it to
155  // a (piecewise-polynomial) high-order mesh.
156  for (int lev = 0; lev < ref_levels; lev++)
157  {
158  mesh->UniformRefinement();
159  }
160  if (mesh->NURBSext)
161  {
162  mesh->SetCurvature(max(order, 1));
163  }
164  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
165 
166  // 5. Define the discontinuous DG finite element space of the given
167  // polynomial order on the refined mesh.
168  DG_FECollection fec(order, dim);
169  FiniteElementSpace fes(mesh, &fec);
170 
171  cout << "Number of unknowns: " << fes.GetVSize() << endl;
172 
173  // 6. Set up and assemble the bilinear and linear forms corresponding to the
174  // DG discretization. The DGTraceIntegrator involves integrals over mesh
175  // interior faces.
179 
180  BilinearForm m(&fes);
182  BilinearForm k(&fes);
183  k.AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
185  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
187  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
188 
189  LinearForm b(&fes);
191  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
192 
193  m.Assemble();
194  m.Finalize();
195  int skip_zeros = 0;
196  k.Assemble(skip_zeros);
197  k.Finalize(skip_zeros);
198  b.Assemble();
199 
200  // 7. Define the initial conditions, save the corresponding grid function to
201  // a file and (optionally) save data in the VisIt format and initialize
202  // GLVis visualization.
203  GridFunction u(&fes);
204  u.ProjectCoefficient(u0);
205 
206  {
207  ofstream omesh("ex9.mesh");
208  omesh.precision(precision);
209  mesh->Print(omesh);
210  ofstream osol("ex9-init.gf");
211  osol.precision(precision);
212  u.Save(osol);
213  }
214 
215  // Create data collection for solution output: either VisItDataCollection for
216  // ascii data files, or SidreDataCollection for binary data files.
217  DataCollection *dc = NULL;
218  if (visit)
219  {
220  if (binary)
221  {
222 #ifdef MFEM_USE_SIDRE
223  dc = new SidreDataCollection("Example9", mesh);
224 #else
225  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
226 #endif
227  }
228  else
229  {
230  dc = new VisItDataCollection("Example9", mesh);
231  dc->SetPrecision(precision);
232  }
233  dc->RegisterField("solution", &u);
234  dc->SetCycle(0);
235  dc->SetTime(0.0);
236  dc->Save();
237  }
238 
239  socketstream sout;
240  if (visualization)
241  {
242  char vishost[] = "localhost";
243  int visport = 19916;
244  sout.open(vishost, visport);
245  if (!sout)
246  {
247  cout << "Unable to connect to GLVis server at "
248  << vishost << ':' << visport << endl;
249  visualization = false;
250  cout << "GLVis visualization disabled.\n";
251  }
252  else
253  {
254  sout.precision(precision);
255  sout << "solution\n" << *mesh << u;
256  sout << "pause\n";
257  sout << flush;
258  cout << "GLVis visualization paused."
259  << " Press space (in the GLVis window) to resume it.\n";
260  }
261  }
262 
263  // 8. Define the time-dependent evolution operator describing the ODE
264  // right-hand side, and perform time-integration (looping over the time
265  // iterations, ti, with a time-step dt).
266  FE_Evolution adv(m.SpMat(), k.SpMat(), b);
267 
268  double t = 0.0;
269  adv.SetTime(t);
270  ode_solver->Init(adv);
271 
272  bool done = false;
273  for (int ti = 0; !done; )
274  {
275  double dt_real = min(dt, t_final - t);
276  ode_solver->Step(u, t, dt_real);
277  ti++;
278 
279  done = (t >= t_final - 1e-8*dt);
280 
281  if (done || ti % vis_steps == 0)
282  {
283  cout << "time step: " << ti << ", time: " << t << endl;
284 
285  if (visualization)
286  {
287  sout << "solution\n" << *mesh << u << flush;
288  }
289 
290  if (visit)
291  {
292  dc->SetCycle(ti);
293  dc->SetTime(t);
294  dc->Save();
295  }
296  }
297  }
298 
299  // 9. Save the final solution. This output can be viewed later using GLVis:
300  // "glvis -m ex9.mesh -g ex9-final.gf".
301  {
302  ofstream osol("ex9-final.gf");
303  osol.precision(precision);
304  u.Save(osol);
305  }
306 
307  // 10. Free the used memory.
308  delete ode_solver;
309  delete dc;
310 
311  return 0;
312 }
313 
314 
315 // Implementation of class FE_Evolution
316 FE_Evolution::FE_Evolution(SparseMatrix &_M, SparseMatrix &_K, const Vector &_b)
317  : TimeDependentOperator(_M.Size()), M(_M), K(_K), b(_b), z(_M.Size())
318 {
319  M_solver.SetPreconditioner(M_prec);
320  M_solver.SetOperator(M);
321 
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);
327 }
328 
329 void FE_Evolution::Mult(const Vector &x, Vector &y) const
330 {
331  // y = M^{-1} (K x + b)
332  K.Mult(x, z);
333  z += b;
334  M_solver.Mult(z, y);
335 }
336 
337 
338 // Velocity coefficient
339 void velocity_function(const Vector &x, Vector &v)
340 {
341  int dim = x.Size();
342 
343  // map to the reference [-1,1] domain
344  Vector X(dim);
345  for (int i = 0; i < dim; i++)
346  {
347  double center = (bb_min[i] + bb_max[i]) * 0.5;
348  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
349  }
350 
351  switch (problem)
352  {
353  case 0:
354  {
355  // Translations in 1D, 2D, and 3D
356  switch (dim)
357  {
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.);
361  break;
362  }
363  break;
364  }
365  case 1:
366  case 2:
367  {
368  // Clockwise rotation in 2D around the origin
369  const double w = M_PI/2;
370  switch (dim)
371  {
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;
375  }
376  break;
377  }
378  case 3:
379  {
380  // Clockwise twisting rotation in 2D around the origin
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.);
383  d = d*d;
384  switch (dim)
385  {
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;
389  }
390  break;
391  }
392  }
393 }
394 
395 // Initial condition
396 double u0_function(const Vector &x)
397 {
398  int dim = x.Size();
399 
400  // map to the reference [-1,1] domain
401  Vector X(dim);
402  for (int i = 0; i < dim; i++)
403  {
404  double center = (bb_min[i] + bb_max[i]) * 0.5;
405  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
406  }
407 
408  switch (problem)
409  {
410  case 0:
411  case 1:
412  {
413  switch (dim)
414  {
415  case 1:
416  return exp(-40.*pow(X(0)-0.5,2));
417  case 2:
418  case 3:
419  {
420  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
421  if (dim == 3)
422  {
423  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
424  rx *= s;
425  ry *= s;
426  }
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;
429  }
430  }
431  }
432  case 2:
433  {
434  double x_ = X(0), y_ = X(1), rho, phi;
435  rho = hypot(x_, y_);
436  phi = atan2(y_, x_);
437  return pow(sin(M_PI*rho),2)*sin(3*phi);
438  }
439  case 3:
440  {
441  const double f = M_PI;
442  return sin(f*X(0))*sin(f*X(1));
443  }
444  }
445  return 0.0;
446 }
447 
448 // Inflow boundary condition (zero for the problems considered in this example)
449 double inflow_function(const Vector &x)
450 {
451  switch (problem)
452  {
453  case 0:
454  case 1:
455  case 2:
456  case 3: return 0.0;
457  }
458  return 0.0;
459 }
Vector bb_max
Definition: ex9.cpp:52
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
int GetVSize() const
Definition: fespace.hpp:163
Conjugate gradient method.
Definition: solvers.hpp:111
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Data type for scaled Jacobi-type smoother of sparse matrix.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.hpp:37
Base abstract class for time dependent operators.
Definition: operator.hpp:140
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Definition: mesh.cpp:88
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:42
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.
Definition: vector.hpp:106
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
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.
Definition: gridfunc.cpp:2185
int dim
Definition: ex3.cpp:47
Data type sparse matrix.
Definition: sparsemat.hpp:38
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6684
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3278
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator.
Definition: linearform.cpp:29
int Dimension() const
Definition: mesh.hpp:611
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Vector bb_min
Definition: ex9.cpp:52
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:147
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:339
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:74
int problem
Definition: ex15.cpp:57
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:134
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:144
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1232
int open(const char hostname[], int port)
double u0_function(const Vector &x)
Definition: ex9.cpp:396
class for C-function coefficient
Vector data type.
Definition: vector.hpp:36
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
The classical forward Euler method.
Definition: ode.hpp:101
double inflow_function(const Vector &x)
Definition: ex9.cpp:449
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:195
virtual void Print(std::ostream &out=std::cout) const
Definition: mesh.hpp:988
bool Good() const
Definition: optparser.hpp:120
alpha (q . grad u, v)