MFEM  v3.0
 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/disc-nurbs.mesh -p 2 -r 3 -dt 0.005 -tf 9
12 // ex9 -m ../data/periodic-square.mesh -p 3 -r 4 -dt 0.0025 -tf 9 -vs 20
13 // ex9 -m ../data/periodic-cube.mesh -p 0 -r 2 -o 2 -dt 0.02 -tf 8
14 //
15 // Description: This example code solves the time-dependent advection equation
16 // du/dt = v.grad(u), where v is a given fluid velocity, and
17 // u0(x)=u(0,x) is a given initial condition.
18 //
19 // The example demonstrates the use of Discontinuous Galerkin (DG)
20 // bilinear forms in MFEM (face integrators), the use of explicit
21 // ODE time integrators, the definition of periodic boundary
22 // conditions through periodic meshes, as well as the use of GLVis
23 // for persistent visualization of a time-evolving solution. The
24 // saving of time-dependent data files for external visualization
25 // with VisIt (visit.llnl.gov) is also illustrated.
26 
27 #include "mfem.hpp"
28 #include <fstream>
29 #include <iostream>
30 #include <algorithm>
31 
32 using namespace std;
33 using namespace mfem;
34 
35 // Choice for the problem setup. The fluid velocity, initial condition and
36 // inflow boundary condition are chosen based on this parameter.
37 int problem;
38 
39 // Velocity coefficient
40 void velocity_function(const Vector &x, Vector &v);
41 
42 // Initial condition
43 double u0_function(Vector &x);
44 
45 // Inflow boundary condition
46 double inflow_function(Vector &x);
47 
48 
54 class FE_Evolution : public TimeDependentOperator
55 {
56 private:
57  SparseMatrix &M, &K;
58  const Vector &b;
59  DSmoother M_prec;
60  CGSolver M_solver;
61 
62  mutable Vector z;
63 
64 public:
65  FE_Evolution(SparseMatrix &_M, SparseMatrix &_K, const Vector &_b);
66 
67  virtual void Mult(const Vector &x, Vector &y) const;
68 
69  virtual ~FE_Evolution() { }
70 };
71 
72 
73 int main(int argc, char *argv[])
74 {
75  // 1. Parse command-line options.
76  problem = 0;
77  const char *mesh_file = "../data/periodic-hexagon.mesh";
78  int ref_levels = 2;
79  int order = 3;
80  int ode_solver_type = 4;
81  double t_final = 10.0;
82  double dt = 0.01;
83  bool visualization = true;
84  bool visit = false;
85  int vis_steps = 5;
86 
87  int precision = 8;
88  cout.precision(precision);
89 
90  OptionsParser args(argc, argv);
91  args.AddOption(&mesh_file, "-m", "--mesh",
92  "Mesh file to use.");
93  args.AddOption(&problem, "-p", "--problem",
94  "Problem setup to use. See options in velocity_function().");
95  args.AddOption(&ref_levels, "-r", "--refine",
96  "Number of times to refine the mesh uniformly.");
97  args.AddOption(&order, "-o", "--order",
98  "Order (degree) of the finite elements.");
99  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
100  "ODE solver: 1 - Forward Euler, 2 - RK2 SSP, 3 - RK3 SSP,"
101  " 4 - RK4, 6 - RK6.");
102  args.AddOption(&t_final, "-tf", "--t-final",
103  "Final time; start time is 0.");
104  args.AddOption(&dt, "-dt", "--time-step",
105  "Time step.");
106  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
107  "--no-visualization",
108  "Enable or disable GLVis visualization.");
109  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
110  "--no-visit-datafiles",
111  "Save data files for VisIt (visit.llnl.gov) visualization.");
112  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
113  "Visualize every n-th timestep.");
114  args.Parse();
115  if (!args.Good())
116  {
117  args.PrintUsage(cout);
118  return 1;
119  }
120  args.PrintOptions(cout);
121 
122  // 2. Read the mesh from the given mesh file. We can handle geometrically
123  // periodic meshes in this code.
124  Mesh *mesh;
125  ifstream imesh(mesh_file);
126  if (!imesh)
127  {
128  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
129  return 2;
130  }
131  mesh = new Mesh(imesh, 1, 1);
132  imesh.close();
133  int dim = mesh->Dimension();
134 
135  // 3. Define the ODE solver used for time integration. Several explicit
136  // Runge-Kutta methods are available.
137  ODESolver *ode_solver = NULL;
138  switch (ode_solver_type)
139  {
140  case 1: ode_solver = new ForwardEulerSolver; break;
141  case 2: ode_solver = new RK2Solver(1.0); break;
142  case 3: ode_solver = new RK3SSPSolver; break;
143  case 4: ode_solver = new RK4Solver; break;
144  case 6: ode_solver = new RK6Solver; break;
145  default:
146  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
147  return 3;
148  }
149 
150  // 4. Refine the mesh to increase the resolution. In this example we do
151  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
152  // command-line parameter. If the mesh is of NURBS type, we convert it to
153  // a (piecewise-polynomial) high-order mesh.
154  for (int lev = 0; lev < ref_levels; lev++)
155  mesh->UniformRefinement();
156 
157  if (mesh->NURBSext)
158  {
159  int mesh_order = std::max(order, 1);
160  FiniteElementCollection *mfec = new H1_FECollection(mesh_order, dim);
161  FiniteElementSpace *mfes = new FiniteElementSpace(mesh, mfec, dim);
162  mesh->SetNodalFESpace(mfes);
163  mesh->GetNodes()->MakeOwner(mfec);
164  }
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  VisItDataCollection visit_dc("Example9", mesh);
216  visit_dc.RegisterField("solution", &u);
217  if (visit)
218  {
219  visit_dc.SetCycle(0);
220  visit_dc.SetTime(0.0);
221  visit_dc.Save();
222  }
223 
224  socketstream sout;
225  if (visualization)
226  {
227  char vishost[] = "localhost";
228  int visport = 19916;
229  sout.open(vishost, visport);
230  if (!sout)
231  {
232  cout << "Unable to connect to GLVis server at "
233  << vishost << ':' << visport << endl;
234  visualization = false;
235  cout << "GLVis visualization disabled.\n";
236  }
237  else
238  {
239  sout.precision(precision);
240  sout << "solution\n" << *mesh << u;
241  sout << "pause\n";
242  sout << flush;
243  cout << "GLVis visualization paused."
244  << " Press space (in the GLVis window) to resume it.\n";
245  }
246  }
247 
248  // 8. Define the time-dependent evolution operator describing the ODE
249  // right-hand side, and perform time-integration (looping over the time
250  // iterations, ti, with a time-step dt).
251  FE_Evolution adv(m.SpMat(), k.SpMat(), b);
252  ode_solver->Init(adv);
253 
254  double t = 0.0;
255  for (int ti = 0; true; )
256  {
257  if (t >= t_final - dt/2)
258  break;
259 
260  ode_solver->Step(u, t, dt);
261  ti++;
262 
263  if (ti % vis_steps == 0)
264  {
265  cout << "time step: " << ti << ", time: " << t << endl;
266 
267  if (visualization)
268  sout << "solution\n" << *mesh << u << flush;
269 
270  if (visit)
271  {
272  visit_dc.SetCycle(ti);
273  visit_dc.SetTime(t);
274  visit_dc.Save();
275  }
276  }
277  }
278 
279  // 9. Save the final solution. This output can be viewed later using GLVis:
280  // "glvis -m ex9.mesh -g ex9-final.gf".
281  {
282  ofstream osol("ex9-final.gf");
283  osol.precision(precision);
284  u.Save(osol);
285  }
286 
287  // 10. Free the used memory.
288  delete ode_solver;
289  delete mesh;
290 
291  return 0;
292 }
293 
294 
295 // Implementation of class FE_Evolution
296 FE_Evolution::FE_Evolution(SparseMatrix &_M, SparseMatrix &_K, const Vector &_b)
297  : TimeDependentOperator(_M.Size()), M(_M), K(_K), b(_b), z(_M.Size())
298 {
299  M_solver.SetPreconditioner(M_prec);
300  M_solver.SetOperator(M);
301 
302  M_solver.iterative_mode = false;
303  M_solver.SetRelTol(1e-9);
304  M_solver.SetAbsTol(0.0);
305  M_solver.SetMaxIter(100);
306  M_solver.SetPrintLevel(0);
307 }
308 
309 void FE_Evolution::Mult(const Vector &x, Vector &y) const
310 {
311  // y = M^{-1} (K x + b)
312  K.Mult(x, z);
313  z += b;
314  M_solver.Mult(z, y);
315 }
316 
317 
318 // Velocity coefficient
319 void velocity_function(const Vector &x, Vector &v)
320 {
321  int dim = x.Size();
322 
323  switch (problem)
324  {
325  case 0:
326  {
327  // Translations in 1D, 2D, and 3D
328  switch (dim)
329  {
330  case 1: v(0) = 1.0; break;
331  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
332  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.); break;
333  }
334  break;
335  }
336  case 1:
337  case 2:
338  {
339  // Clockwise rotation in 2D around the origin
340  const double w = M_PI/2;
341  switch (dim)
342  {
343  case 1: v(0) = 1.0; break;
344  case 2: v(0) = w*x(1); v(1) = -w*x(0); break;
345  case 3: v(0) = w*x(1); v(1) = -w*x(0); v(2) = 0.0; break;
346  }
347  break;
348  }
349  case 3:
350  {
351  // Clockwise twisting rotation in 2D around the origin
352  const double w = M_PI/2;
353  double d = max((x(0)+1.)*(1.-x(0)),0.) * max((x(1)+1.)*(1.-x(1)),0.);
354  d = d*d;
355  switch (dim)
356  {
357  case 1: v(0) = 1.0; break;
358  case 2: v(0) = d*w*x(1); v(1) = -d*w*x(0); break;
359  case 3: v(0) = d*w*x(1); v(1) = -d*w*x(0); v(2) = 0.0; break;
360  }
361  break;
362  }
363  }
364 }
365 
366 // Initial condition
367 double u0_function(Vector &x)
368 {
369  int dim = x.Size();
370 
371  switch (problem)
372  {
373  case 0:
374  case 1:
375  {
376  switch (dim)
377  {
378  case 1:
379  return exp(-40.*pow(x(0)-0.5,2));
380  case 2:
381  case 3:
382  {
383  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
384  if (dim == 3)
385  {
386  const double s = (1. + 0.25*cos(2*M_PI*x(2)));
387  rx *= s;
388  ry *= s;
389  }
390  return ( erfc(w*(x(0)-cx-rx))*erfc(-w*(x(0)-cx+rx)) *
391  erfc(w*(x(1)-cy-ry))*erfc(-w*(x(1)-cy+ry)) )/16;
392  }
393  }
394  }
395  case 2:
396  {
397  const double r = sqrt(8.);
398  double x_ = x(0), y_ = x(1), rho, phi;
399  rho = hypot(x_, y_) / r;
400  phi = atan2(y_, x_);
401  return pow(sin(M_PI*rho),2)*sin(3*phi);
402  }
403  case 3:
404  {
405  const double f = M_PI;
406  return sin(f*x(0))*sin(f*x(1));
407  }
408  }
409  return 0.0;
410 }
411 
412 // Inflow boundary condition (zero for the problems considered in this example)
414 {
415  switch (problem)
416  {
417  case 0:
418  case 1:
419  case 2:
420  case 3: return 0.0;
421  }
422  return 0.0;
423 }
int GetVSize() const
Definition: fespace.hpp:151
Conjugate gradient method.
Definition: solvers.hpp:109
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:26
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 RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Base abstract class for time dependent operators: (x,t) -&gt; f(x,t)
Definition: operator.hpp:68
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:351
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
virtual void Step(Vector &x, double &t, double &dt)=0
int Size() const
Returns the size of the vector.
Definition: vector.hpp:76
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
double inflow_function(Vector &x)
Definition: ex9.cpp:413
virtual void Save()
Save the collection and a VisIt root file.
virtual void Init(TimeDependentOperator &_f)
Definition: ode.hpp:30
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:1797
Data type sparse matrix.
Definition: sparsemat.hpp:38
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
Data collection with VisIt I/O routines.
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator.
Definition: linearform.cpp:29
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3066
int problem
Definition: ex9.cpp:37
int Dimension() const
Definition: mesh.hpp:417
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:385
void SetTime(double t)
Set physical time (for time-dependent simulations)
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Definition: mesh.cpp:7136
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:85
int main(int argc, char *argv[])
Definition: ex1.cpp:39
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:319
Abstract finite element space.
Definition: fespace.hpp:61
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
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:72
double u0_function(Vector &x)
Definition: ex9.cpp:367
NURBSExtension * NURBSext
Definition: mesh.hpp:307
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:266
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:967
int open(const char hostname[], int port)
class for C-function coefficient
Vector data type.
Definition: vector.hpp:29
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:4948
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:52
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
const SparseMatrix & SpMat() const
Returns a reference to the sparse martix.
The classical forward Euler method.
Definition: ode.hpp:39
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:83
bool Good() const
Definition: optparser.hpp:120
alpha (q . grad u, v)
Definition: bilininteg.hpp:239