MFEM  v4.1.0
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 // Nonlinear Constrained Optimization Modification
3 //
4 // Compile with: make ex9
5 //
6 // Sample runs:
7 // ex9 -m ../../data/periodic-segment.mesh -r 3 -p 0 -o 2 -dt 0.002 -opt 1
8 // ex9 -m ../../data/periodic-segment.mesh -r 3 -p 0 -o 2 -dt 0.002 -opt 2
9 //
10 // ex9 -m ../../data/periodic-square.mesh -p 0 -r 2 -dt 0.01 -tf 10 -opt 1
11 // ex9 -m ../../data/periodic-square.mesh -p 0 -r 2 -dt 0.01 -tf 10 -opt 2
12 //
13 // ex9 -m ../../data/periodic-square.mesh -p 1 -r 2 -dt 0.005 -tf 9 -opt 1
14 // ex9 -m ../../data/periodic-square.mesh -p 1 -r 2 -dt 0.005 -tf 9 -opt 2
15 //
16 // ex9 -m ../../data/amr-quad.mesh -p 1 -r 1 -dt 0.002 -tf 9 -opt 1
17 // ex9 -m ../../data/amr-quad.mesh -p 1 -r 1 -dt 0.002 -tf 9 -opt 2
18 //
19 // ex9 -m ../../data/disc-nurbs.mesh -p 1 -r 2 -dt 0.005 -tf 9 -opt 1
20 // ex9 -m ../../data/disc-nurbs.mesh -p 1 -r 2 -dt 0.005 -tf 9 -opt 2
21 //
22 // ex9 -m ../../data/disc-nurbs.mesh -p 2 -r 2 -dt 0.01 -tf 9 -opt 1
23 // ex9 -m ../../data/disc-nurbs.mesh -p 2 -r 2 -dt 0.01 -tf 9 -opt 2
24 //
25 // ex9 -m ../../data/periodic-square.mesh -p 3 -r 3 -dt 0.0025 -tf 9 -opt 1
26 // ex9 -m ../../data/periodic-square.mesh -p 3 -r 3 -dt 0.0025 -tf 9 -opt 2
27 //
28 // ex9 -m ../../data/periodic-cube.mesh -p 0 -r 2 -o 2 -dt 0.02 -tf 8 -opt 1
29 // ex9 -m ../../data/periodic-cube.mesh -p 0 -r 2 -o 2 -dt 0.02 -tf 8 -opt 2
30 
31 // Description: This example modifies the standard MFEM ex9 by adding nonlinear
32 // constrained optimization capabilities through the SLBQP and
33 // HIOP solvers. It demonstrates how a user can define a custom
34 // class OptimizationProblem that includes linear/nonlinear
35 // equality/inequality constraints. This optimization is applied
36 // as post-processing to the solution of the transport equation.
37 //
38 // Description of ex9:
39 // This example code solves the time-dependent advection equation
40 // du/dt + v.grad(u) = 0, where v is a given fluid velocity, and
41 // u0(x)=u(0,x) is a given initial condition.
42 //
43 // The example demonstrates the use of Discontinuous Galerkin (DG)
44 // bilinear forms in MFEM (face integrators), the use of explicit
45 // ODE time integrators, the definition of periodic boundary
46 // conditions through periodic meshes, as well as the use of GLVis
47 // for persistent visualization of a time-evolving solution. The
48 // saving of time-dependent data files for external visualization
49 // with VisIt (visit.llnl.gov) is also illustrated.
50 
51 #include "mfem.hpp"
52 #include <fstream>
53 #include <iostream>
54 
55 using namespace std;
56 using namespace mfem;
57 
58 // Choice for the problem setup. The fluid velocity, initial condition and
59 // inflow boundary condition are chosen based on this parameter.
60 int problem;
61 
62 // Nonlinear optimizer.
64 
65 // Velocity coefficient
66 bool invert_velocity = false;
67 void velocity_function(const Vector &x, Vector &v);
68 
69 // Initial condition
70 double u0_function(const Vector &x);
71 
72 // Inflow boundary condition
73 double inflow_function(const Vector &x);
74 
75 // Mesh bounding box
77 
78 /// Computes C(x) = sum w_i x_i, where w is a given Vector.
79 class LinearScaleOperator : public Operator
80 {
81 private:
82  const Vector &w;
83  mutable DenseMatrix grad;
84 
85 public:
86  LinearScaleOperator(const Vector &weight)
87  : Operator(1, weight.Size()), w(weight), grad(1, width)
88  {
89  for (int i = 0; i < width; i++) { grad(0, i) = w(i); }
90  }
91 
92  virtual void Mult(const Vector &x, Vector &y) const
93  {
94  y(0) = w * x;
95  }
96 
97  virtual Operator &GetGradient(const Vector &x) const
98  {
99  return grad;
100  }
101 };
102 
103 /// Nonlinear monotone bounded operator to test nonlinear ineq constraints.
104 /// Computes D(x) = tanh(sum(x_i)).
105 class TanhSumOperator : public Operator
106 {
107 private:
108  mutable DenseMatrix grad;
109 
110 public:
111  TanhSumOperator(int size) : Operator(1, size), grad(1, width) { }
112 
113  virtual void Mult(const Vector &x, Vector &y) const
114  {
115  y(0) = std::tanh(x.Sum());
116  }
117 
118  virtual Operator &GetGradient(const Vector &x) const
119  {
120  const double ts = std::tanh(x.Sum());
121  const double dtanh = 1.0 - ts * ts;
122  for (int i = 0; i < width; i++) { grad(0, i) = dtanh; }
123  return grad;
124  }
125 };
126 
127 /** Monotone and conservative a-posteriori correction for transport solutions:
128  * Find x that minimizes 0.5 || x - x_HO ||^2, subject to
129  * sum w_i x_i = mass,
130  * tanh(sum(x_i_min)) <= tanh(sum(x_i)) <= tanh(sum(x_i_max)),
131  * x_i_min <= x_i <= x_i_max,
132  */
133 class OptimizedTransportProblem : public OptimizationProblem
134 {
135 private:
136  const Vector &x_HO;
137  Vector massvec, d_lo, d_hi;
138  const LinearScaleOperator LSoper;
139  const TanhSumOperator TSoper;
140 
141 public:
142  OptimizedTransportProblem(const Vector &xho, const Vector &w, double mass,
143  const Vector &xmin, const Vector &xmax)
144  : OptimizationProblem(xho.Size(), NULL, NULL),
145  x_HO(xho), massvec(1), d_lo(1), d_hi(1),
146  LSoper(w), TSoper(w.Size())
147  {
148  C = &LSoper;
149  massvec(0) = mass;
150  SetEqualityConstraint(massvec);
151 
152  D = &TSoper;
153  d_lo(0) = std::tanh(xmin.Sum());
154  d_hi(0) = std::tanh(xmax.Sum());
155  MFEM_ASSERT(d_lo(0) < d_hi(0),
156  "The bounds produce an infeasible optimization problem");
157  SetInequalityConstraint(d_lo, d_hi);
158 
159  SetSolutionBounds(xmin, xmax);
160  }
161 
162  virtual double CalcObjective(const Vector &x) const
163  {
164  double res = 0.0;
165  for (int i = 0; i < input_size; i++)
166  {
167  const double d = x(i) - x_HO(i);
168  res += d * d;
169  }
170  return 0.5 * res;
171  }
172 
173  virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
174  {
175  for (int i = 0; i < input_size; i++) { grad(i) = x(i) - x_HO(i); }
176  }
177 };
178 
179 
180 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
181  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
182  and advection matrices, and b describes the flow on the boundary. This can
183  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
184  used to evaluate the right-hand side. */
185 class FE_Evolution : public TimeDependentOperator
186 {
187 private:
188  SparseMatrix &M, &K;
189  const Vector &b;
190  DSmoother M_prec;
191  CGSolver M_solver;
192 
193  mutable Vector z;
194 
195  double dt;
196  BilinearForm &bf;
197  Vector &M_rowsums;
198 
199 public:
200  FE_Evolution(SparseMatrix &_M, SparseMatrix &_K, const Vector &_b,
201  BilinearForm &_bf, Vector &M_rs);
202 
203  void SetTimeStep(double _dt) { dt = _dt; }
204  void SetK(SparseMatrix &_K) { K = _K; }
205  virtual void Mult(const Vector &x, Vector &y) const;
206 
207  virtual ~FE_Evolution() { }
208 };
209 
210 
211 int main(int argc, char *argv[])
212 {
213  // 1. Parse command-line options.
214  problem = 0;
215  optimizer_type = 1;
216  const char *mesh_file = "../../data/periodic-hexagon.mesh";
217  int ref_levels = 2;
218  int order = 3;
219  int ode_solver_type = 3;
220  double t_final = 1.0;
221  double dt = 0.01;
222  bool visualization = true;
223  bool visit = false;
224  bool binary = false;
225  int vis_steps = 5;
226 
227  int precision = 8;
228  cout.precision(precision);
229 
230  OptionsParser args(argc, argv);
231  args.AddOption(&mesh_file, "-m", "--mesh",
232  "Mesh file to use.");
233  args.AddOption(&problem, "-p", "--problem",
234  "Problem setup to use. See options in velocity_function().");
235  args.AddOption(&ref_levels, "-r", "--refine",
236  "Number of times to refine the mesh uniformly.");
237  args.AddOption(&order, "-o", "--order",
238  "Order (degree) of the finite elements.");
239  args.AddOption(&optimizer_type, "-opt", "--optimizer",
240  "Nonlinear optimizer: 1 - SLBQP,\n\t"
241  " 2 - HIOP.");
242  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
243  "ODE solver: 1 - Forward Euler,\n\t"
244  " 2 - RK2 SSP, 3 - RK3 SSP, 4 - RK4, 6 - RK6.");
245  args.AddOption(&t_final, "-tf", "--t-final",
246  "Final time; start time is 0.");
247  args.AddOption(&dt, "-dt", "--time-step",
248  "Time step.");
249  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
250  "--no-visualization",
251  "Enable or disable GLVis visualization.");
252  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
253  "--no-visit-datafiles",
254  "Save data files for VisIt (visit.llnl.gov) visualization.");
255  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
256  "--ascii-datafiles",
257  "Use binary (Sidre) or ascii format for VisIt data files.");
258  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
259  "Visualize every n-th timestep.");
260  args.Parse();
261  if (!args.Good())
262  {
263  args.PrintUsage(cout);
264  return 1;
265  }
266  args.PrintOptions(cout);
267 
268  // 2. Read the mesh from the given mesh file. We can handle geometrically
269  // periodic meshes in this code.
270  Mesh *mesh = new Mesh(mesh_file, 1, 1);
271  int dim = mesh->Dimension();
272 
273  // 3. Define the ODE solver used for time integration. Several explicit
274  // Runge-Kutta methods are available.
275  ODESolver *ode_solver = NULL;
276  switch (ode_solver_type)
277  {
278  case 1: ode_solver = new ForwardEulerSolver; break;
279  case 2: ode_solver = new RK2Solver(1.0); break;
280  case 3: ode_solver = new RK3SSPSolver; break;
281  case 4: ode_solver = new RK4Solver; break;
282  case 6: ode_solver = new RK6Solver; break;
283  default:
284  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
285  delete mesh;
286  return 3;
287  }
288 
289  // 4. Refine the mesh to increase the resolution. In this example we do
290  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
291  // command-line parameter. If the mesh is of NURBS type, we convert it to
292  // a (piecewise-polynomial) high-order mesh.
293  for (int lev = 0; lev < ref_levels; lev++)
294  {
295  mesh->UniformRefinement();
296  }
297  if (mesh->NURBSext)
298  {
299  mesh->SetCurvature(max(order, 1));
300  }
301  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
302 
303  // 5. Define the discontinuous DG finite element space of the given
304  // polynomial order on the refined mesh.
305  DG_FECollection fec(order, dim, BasisType::Positive);
306  FiniteElementSpace fes(mesh, &fec);
307 
308  cout << "Number of unknowns: " << fes.GetVSize() << endl;
309 
310  // 6. Set up and assemble the bilinear and linear forms corresponding to the
311  // DG discretization. The DGTraceIntegrator involves integrals over mesh
312  // interior faces.
316 
317  BilinearForm m(&fes);
319  BilinearForm k(&fes);
320  k.AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
322  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
324  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
325 
326  LinearForm b(&fes);
328  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
329 
330  m.Assemble();
331  m.Finalize();
332  int skip_zeros = 0;
333  k.Assemble(skip_zeros);
334  k.Finalize(skip_zeros);
335  b.Assemble();
336 
337  // 7. Define the initial conditions, save the corresponding grid function to
338  // a file and (optionally) save data in the VisIt format and initialize
339  // GLVis visualization.
340  GridFunction u(&fes);
341  u.ProjectCoefficient(u0);
342 
343  {
344  ofstream omesh("ex9.mesh");
345  omesh.precision(precision);
346  mesh->Print(omesh);
347  ofstream osol("ex9-init.gf");
348  osol.precision(precision);
349  u.Save(osol);
350  }
351 
352  // Create data collection for solution output: either VisItDataCollection for
353  // ascii data files, or SidreDataCollection for binary data files.
354  DataCollection *dc = NULL;
355  if (visit)
356  {
357  if (binary)
358  {
359 #ifdef MFEM_USE_SIDRE
360  dc = new SidreDataCollection("Example9", mesh);
361 #else
362  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
363 #endif
364  }
365  else
366  {
367  dc = new VisItDataCollection("Example9", mesh);
368  dc->SetPrecision(precision);
369  }
370  dc->RegisterField("solution", &u);
371  dc->SetCycle(0);
372  dc->SetTime(0.0);
373  dc->Save();
374  }
375 
376  socketstream sout;
377  if (visualization)
378  {
379  char vishost[] = "localhost";
380  int visport = 19916;
381  sout.open(vishost, visport);
382  if (!sout)
383  {
384  cout << "Unable to connect to GLVis server at "
385  << vishost << ':' << visport << endl;
386  visualization = false;
387  cout << "GLVis visualization disabled.\n";
388  }
389  else
390  {
391  sout.precision(precision);
392  sout << "solution\n" << *mesh << u;
393  sout << "pause\n";
394  sout << flush;
395  cout << "GLVis visualization paused."
396  << " Press space (in the GLVis window) to resume it.\n";
397  }
398  }
399 
400  Vector M_rowsums(m.Size());
401  m.SpMat().GetRowSums(M_rowsums);
402 
403  // 8. Define the time-dependent evolution operator describing the ODE
404  // right-hand side, and perform time-integration (looping over the time
405  // iterations, ti, with a time-step dt).
406  FE_Evolution adv(m.SpMat(), k.SpMat(), b, k, M_rowsums);
407 
408  double t = 0.0;
409  adv.SetTime(t);
410  ode_solver->Init(adv);
411 
412  // Compute initial volume.
413  const double vol0 = M_rowsums * u;
414 
415  bool done = false;
416  for (int ti = 0; !done; )
417  {
418  double dt_real = min(dt, t_final - t);
419  adv.SetTimeStep(dt_real);
420  ode_solver->Step(u, t, dt_real);
421  ti++;
422 
423  done = (t >= t_final - 1e-8*dt);
424 
425  if (done || ti % vis_steps == 0)
426  {
427  cout << "time step: " << ti << ", time: " << t << endl;
428 
429  if (visualization)
430  {
431  sout << "solution\n" << *mesh << u << flush;
432  }
433 
434  if (visit)
435  {
436  dc->SetCycle(ti);
437  dc->SetTime(t);
438  dc->Save();
439  }
440  }
441  }
442 
443  // Print the error vs exact solution.
444  const double max_error = u.ComputeMaxError(u0),
445  l1_error = u.ComputeL1Error(u0),
446  l2_error = u.ComputeL2Error(u0);
447  std::cout << "Linf error = " << max_error << endl
448  << "L1 error = " << l1_error << endl
449  << "L2 error = " << l2_error << endl;
450 
451  // Print error in volume.
452  const double vol = M_rowsums * u;
453  std::cout << "Vol error = " << vol - vol0 << endl;
454 
455  // 9. Save the final solution. This output can be viewed later using GLVis:
456  // "glvis -m ex9.mesh -g ex9-final.gf".
457  {
458  ofstream osol("ex9-final.gf");
459  osol.precision(precision);
460  u.Save(osol);
461  }
462 
463  // 10. Free the used memory.
464  delete ode_solver;
465  delete dc;
466  delete mesh;
467 
468  return 0;
469 }
470 
471 
472 // Implementation of class FE_Evolution
474  const Vector &_b, BilinearForm &_bf, Vector &M_rs)
475  : TimeDependentOperator(_M.Size()),
476  M(_M), K(_K), b(_b), M_prec(), M_solver(), z(_M.Size()),
477  bf(_bf), M_rowsums(M_rs)
478 {
479  M_solver.SetPreconditioner(M_prec);
480  M_solver.SetOperator(M);
481 
482  M_solver.iterative_mode = false;
483  M_solver.SetRelTol(1e-9);
484  M_solver.SetAbsTol(0.0);
485  M_solver.SetMaxIter(100);
486  M_solver.SetPrintLevel(0);
487 }
488 
489 void FE_Evolution::Mult(const Vector &x, Vector &y) const
490 {
491  // Compute bounds y_min, y_max for y from x on the ldofs.
492  const int dofs = x.Size();
493  Vector y_min(dofs), y_max(dofs);
494  const int *In = bf.SpMat().GetI(), *Jn = bf.SpMat().GetJ();
495  for (int i = 0, k = 0; i < dofs; i++)
496  {
497  double x_i_min = +std::numeric_limits<double>::infinity();
498  double x_i_max = -std::numeric_limits<double>::infinity();
499  for (int end = In[i+1]; k < end; k++)
500  {
501  const int j = Jn[k];
502  if (x(j) > x_i_max) { x_i_max = x(j); }
503  if (x(j) < x_i_min) { x_i_min = x(j); }
504  }
505  y_min(i) = x_i_min;
506  y_max(i) = x_i_max;
507  }
508  for (int i = 0; i < dofs; i++)
509  {
510  y_min(i) = (y_min(i) - x(i) ) / dt;
511  y_max(i) = (y_max(i) - x(i) ) / dt;
512  }
513 
514  // Compute the high-order solution y = M^{-1} (K x + b).
515  K.Mult(x, z);
516  z += b;
517  M_solver.Mult(z, y);
518 
519  // The solution y is an increment; it should not introduce new mass.
520  const double mass_y = 0.0;
521 
522  // Perform optimization.
523  Vector y_out(dofs);
524  const int max_iter = 500;
525  const double rtol = 1.e-7;
526  double atol = 1.e-7;
527 
528  OptimizationSolver *optsolver = NULL;
529  if (optimizer_type == 2)
530  {
531 #ifdef MFEM_USE_HIOP
532  HiopNlpOptimizer *tmp_opt_ptr = new HiopNlpOptimizer();
533  optsolver = tmp_opt_ptr;
534 #else
535  MFEM_ABORT("MFEM is not built with HiOp support!");
536 #endif
537  }
538  else
539  {
540  SLBQPOptimizer *slbqp = new SLBQPOptimizer();
541  slbqp->SetBounds(y_min, y_max);
542  slbqp->SetLinearConstraint(M_rowsums, mass_y);
543  atol = 1.e-15;
544  optsolver = slbqp;
545  }
546 
547  OptimizedTransportProblem ot_prob(y, M_rowsums, mass_y, y_min, y_max);
548  optsolver->SetOptimizationProblem(ot_prob);
549 
550  optsolver->SetMaxIter(max_iter);
551  optsolver->SetAbsTol(atol);
552  optsolver->SetRelTol(rtol);
553  optsolver->SetPrintLevel(0);
554  optsolver->Mult(y, y_out);
555 
556  y = y_out;
557 
558  delete optsolver;
559 }
560 
561 
562 // Velocity coefficient
563 void velocity_function(const Vector &x, Vector &v)
564 {
565  int dim = x.Size();
566 
567  // map to the reference [-1,1] domain
568  Vector X(dim);
569  for (int i = 0; i < dim; i++)
570  {
571  double center = (bb_min[i] + bb_max[i]) * 0.5;
572  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
573  }
574 
575  switch (problem)
576  {
577  case 0:
578  {
579  // Translations in 1D, 2D, and 3D
580  switch (dim)
581  {
582  case 1: v(0) = (invert_velocity) ? -1.0 : 1.0; break;
583  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
584  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
585  break;
586  }
587  break;
588  }
589  case 1:
590  case 2:
591  {
592  // Clockwise rotation in 2D around the origin
593  const double w = M_PI/2;
594  switch (dim)
595  {
596  case 1: v(0) = 1.0; break;
597  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
598  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
599  }
600  break;
601  }
602  case 3:
603  {
604  // Clockwise twisting rotation in 2D around the origin
605  const double w = M_PI/2;
606  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
607  d = d*d;
608  switch (dim)
609  {
610  case 1: v(0) = 1.0; break;
611  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
612  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
613  }
614  break;
615  }
616  }
617 }
618 
619 // Initial condition
620 double u0_function(const Vector &x)
621 {
622  int dim = x.Size();
623 
624  // map to the reference [-1,1] domain
625  Vector X(dim);
626  for (int i = 0; i < dim; i++)
627  {
628  double center = (bb_min[i] + bb_max[i]) * 0.5;
629  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
630  }
631 
632  switch (problem)
633  {
634  case 0:
635  case 1:
636  {
637  switch (dim)
638  {
639  case 1:
640  return (X(0) > -0.15 && X(0) < 0.15) ? 1.0 : 0.0;
641  //return exp(-40.*pow(X(0)-0.0,2));
642  case 2:
643  case 3:
644  {
645  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
646  if (dim == 3)
647  {
648  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
649  rx *= s;
650  ry *= s;
651  }
652  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
653  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
654  }
655  }
656  }
657  case 2:
658  {
659  double x_ = X(0), y_ = X(1), rho, phi;
660  rho = hypot(x_, y_);
661  phi = atan2(y_, x_);
662  return pow(sin(M_PI*rho),2)*sin(3*phi);
663  }
664  case 3:
665  {
666  const double f = M_PI;
667  return sin(f*X(0))*sin(f*X(1));
668  }
669  }
670  return 0.0;
671 }
672 
673 // Inflow boundary condition (zero for the problems considered in this example)
674 double inflow_function(const Vector &x)
675 {
676  switch (problem)
677  {
678  case 0:
679  case 1:
680  case 2:
681  case 3: return 0.0;
682  }
683  return 0.0;
684 }
Vector bb_max
Definition: ex9.cpp:60
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:381
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1188
Conjugate gradient method.
Definition: solvers.hpp:152
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
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 Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:358
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
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:120
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:146
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
virtual void SetTime(const double _t)
Set the current time.
Definition: operator.hpp:303
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
int * GetI()
Return the array I.
Definition: sparsemat.hpp:141
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:504
int main(int argc, char *argv[])
Definition: ex1.cpp:62
FE_Evolution(FiniteElementSpace &_vfes, Operator &_A, SparseMatrix &_Aflux)
Definition: ex18.hpp:104
Data collection with Sidre routines following the Conduit mesh blueprint specification.
int optimizer_type
Definition: ex9.cpp:63
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2731
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetLinearConstraint(const Vector &_w, double _a)
Definition: solvers.cpp:1632
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
Data type sparse matrix.
Definition: sparsemat.hpp:40
virtual void Save()
Save the collection to disk.
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
int Size() const
Get the size of the BilinearForm as a square matrix.
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3924
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:66
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)
void SetTimeStep(double _dt)
Definition: ex9.cpp:203
Vector bb_min
Definition: ex9.cpp:60
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:145
void SetAbsTol(double atol)
Definition: solvers.hpp:64
void SetRelTol(double rtol)
Definition: solvers.hpp:63
void velocity_function(const Vector &x, Vector &v)
Definition: ex9.cpp:478
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
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
int problem
Definition: ex15.cpp:54
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:132
virtual ~FE_Evolution()
Definition: ex9.cpp:207
Adapts the HiOp functionality to the MFEM OptimizationSolver interface.
Definition: hiop.hpp:172
virtual void Mult(const Vector &x, Vector &y) const
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x...
Definition: ex18.hpp:130
void SetBounds(const Vector &_lo, const Vector &_hi)
Definition: solvers.cpp:1626
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:191
Abstract solver for OptimizationProblems.
Definition: solvers.hpp:394
bool invert_velocity
Definition: ex9.cpp:66
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual void Mult(const Vector &xt, Vector &x) const =0
Operator application: y=A(x).
int dim
Definition: ex24.cpp:43
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1685
int open(const char hostname[], int port)
double u0_function(const Vector &x)
Definition: ex9.cpp:535
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:166
class for C-function coefficient
void SetK(SparseMatrix &_K)
Definition: ex9.cpp:204
Vector data type.
Definition: vector.hpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
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 matrix.
The classical forward Euler method.
Definition: ode.hpp:99
Abstract operator.
Definition: operator.hpp:24
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
Definition: solvers.hpp:408
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:816
double inflow_function(const Vector &x)
Definition: ex9.cpp:588
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
void GetRowSums(Vector &x) const
For all i compute .
Definition: sparsemat.cpp:823
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:143
bool Good() const
Definition: optparser.hpp:125
alpha (q . grad u, v)