MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex9p.cpp
Go to the documentation of this file.
1 // MFEM Example 9 - Parallel Version
2 // Nonlinear Constrained Optimization Modification
3 //
4 // Compile with: make ex9p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex9p -m ../../data/periodic-segment.mesh -rs 3 -p 0 -o 2 -dt 0.002 -opt 1
8 // mpirun -np 4 ex9p -m ../../data/periodic-segment.mesh -rs 3 -p 0 -o 2 -dt 0.002 -opt 2
9 //
10 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 0 -rs 2 -dt 0.01 -tf 10 -opt 1
11 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 0 -rs 2 -dt 0.01 -tf 10 -opt 2
12 //
13 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 1 -rs 2 -dt 0.005 -tf 9 -opt 1
14 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 1 -rs 2 -dt 0.005 -tf 9 -opt 2
15 //
16 // mpirun -np 4 ex9p -m ../../data/amr-quad.mesh -p 1 -rs 1 -dt 0.002 -tf 9 -opt 1
17 // mpirun -np 4 ex9p -m ../../data/amr-quad.mesh -p 1 -rs 1 -dt 0.002 -tf 9 -opt 2
18 //
19 // mpirun -np 4 ex9p -m ../../data/disc-nurbs.mesh -p 1 -rs 2 -dt 0.005 -tf 9 -opt 1
20 // mpirun -np 4 ex9p -m ../../data/disc-nurbs.mesh -p 1 -rs 2 -dt 0.005 -tf 9 -opt 2
21 //
22 // mpirun -np 4 ex9p -m ../../data/disc-nurbs.mesh -p 2 -rs 2 -dt 0.01 -tf 9 -opt 1
23 // mpirun -np 4 ex9p -m ../../data/disc-nurbs.mesh -p 2 -rs 2 -dt 0.01 -tf 9 -opt 2
24 //
25 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 3 -rs 3 -dt 0.0025 -tf 9 -opt 1
26 // mpirun -np 4 ex9p -m ../../data/periodic-square.mesh -p 3 -rs 3 -dt 0.0025 -tf 9 -opt 2
27 //
28 // mpirun -np 4 ex9p -m ../../data/periodic-cube.mesh -p 0 -rs 2 -o 2 -dt 0.02 -tf 8 -opt 1
29 // mpirun -np 4 ex9p -m ../../data/periodic-cube.mesh -p 0 -rs 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:
83  // Local weights.
84  const Vector &w;
85  // Gradient for the tdofs.
86  mutable DenseMatrix grad;
87 
88 public:
89  LinearScaleOperator(ParFiniteElementSpace &space, const Vector &weight)
90  : Operator(1, space.TrueVSize()),
91  pfes(space), w(weight), grad(1, width)
92  {
93  Vector w_glob(width);
94  pfes.Dof_TrueDof_Matrix()->MultTranspose(w, w_glob);
95  for (int i = 0; i < width; i++) { grad(0, i) = w_glob(i); }
96  }
97 
98  virtual void Mult(const Vector &x, Vector &y) const
99  {
100  Vector x_loc(w.Size());
101  pfes.GetProlongationMatrix()->Mult(x, x_loc);
102  const double loc_res = w * x_loc;
103  MPI_Allreduce(&loc_res, &y(0), 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
104  }
105 
106  virtual Operator &GetGradient(const Vector &x) const
107  {
108  return grad;
109  }
110 };
111 
112 /// Nonlinear monotone bounded operator to test nonlinear ineq constraints.
113 /// Computes D(x) = tanh(sum(x_i)).
114 class TanhSumOperator : public Operator
115 {
116 private:
117  // Gradient for the tdofs.
118  mutable DenseMatrix grad;
119 
120 public:
121  TanhSumOperator(ParFiniteElementSpace &space)
122  : Operator(1, space.TrueVSize()), grad(1, width) { }
123 
124  virtual void Mult(const Vector &x, Vector &y) const
125  {
126  double sum_loc = x.Sum();
127  MPI_Allreduce(&sum_loc, &y(0), 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
128  y(0) = std::tanh(y(0));
129  }
130 
131  virtual Operator &GetGradient(const Vector &x) const
132  {
133  double sum_loc = x.Sum();
134  double dtanh;
135  MPI_Allreduce(&sum_loc, &dtanh, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
136  dtanh = 1.0 - pow(std::tanh(dtanh), 2);
137 
138  for (int i = 0; i < width; i++) { grad(0, i) = dtanh; }
139  return grad;
140  }
141 };
142 
143 
144 /** Monotone and conservative a-posteriori correction for transport solutions:
145  * Find x that minimizes 0.5 || x - x_HO ||^2, subject to
146  * sum w_i x_i = mass,
147  * tanh(sum(x_i_min)) <= tanh(sum(x_i)) <= tanh(sum(x_i_max)),
148  * x_i_min <= x_i <= x_i_max,
149  */
150 class OptimizedTransportProblem : public OptimizationProblem
151 {
152 private:
153  const Vector &x_HO;
154  Vector massvec, d_lo, d_hi;
155  const LinearScaleOperator LSoper;
156  const TanhSumOperator TSoper;
157 
158 public:
159  OptimizedTransportProblem(ParFiniteElementSpace &space,
160  const Vector &xho, const Vector &w, double mass,
161  const Vector &xmin, const Vector &xmax)
162  : OptimizationProblem(xho.Size(), NULL, NULL),
163  x_HO(xho), massvec(1), d_lo(1), d_hi(1),
164  LSoper(space, w), TSoper(space)
165  {
166  C = &LSoper;
167  massvec(0) = mass;
168  SetEqualityConstraint(massvec);
169 
170  D = &TSoper;
171  double lsums[2], gsums[2];
172  lsums[0] = xmin.Sum();
173  lsums[1] = xmax.Sum();
174  MPI_Allreduce(lsums, gsums, 2, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
175  d_lo(0) = std::tanh(gsums[0]);
176  d_hi(0) = std::tanh(gsums[1]);
177  MFEM_ASSERT(d_lo(0) < d_hi(0),
178  "The bounds produce an infeasible optimization problem");
179  SetInequalityConstraint(d_lo, d_hi);
180 
181  SetSolutionBounds(xmin, xmax);
182  }
183 
184  virtual double CalcObjective(const Vector &x) const
185  {
186  double loc_res = 0.0;
187  for (int i = 0; i < input_size; i++)
188  {
189  const double d = x(i) - x_HO(i);
190  loc_res += d * d;
191  }
192  loc_res *= 0.5;
193  double res;
194  MPI_Allreduce(&loc_res, &res, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
195  return res;
196  }
197 
198  virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
199  {
200  for (int i = 0; i < input_size; i++) { grad(i) = x(i) - x_HO(i); }
201  }
202 };
203 
204 
205 /** A time-dependent operator for the right-hand side of the ODE. The DG weak
206  form of du/dt = -v.grad(u) is M du/dt = K u + b, where M and K are the mass
207  and advection matrices, and b describes the flow on the boundary. This can
208  be written as a general ODE, du/dt = M^{-1} (K u + b), and this class is
209  used to evaluate the right-hand side. */
210 class FE_Evolution : public TimeDependentOperator
211 {
212 private:
213  HypreParMatrix &M, &K;
214  const Vector &b;
215  HypreSmoother M_prec;
216  CGSolver M_solver;
217 
218  mutable Vector z;
219 
220  double dt;
221  ParBilinearForm &pbf;
222  Vector &M_rowsums;
223 
224 public:
226  const Vector &_b, ParBilinearForm &_pbf, Vector &M_rs);
227 
228  void SetTimeStep(double _dt) { dt = _dt; }
229  void SetK(HypreParMatrix &_K) { K = _K; }
230  virtual void Mult(const Vector &x, Vector &y) const;
231 
232  virtual ~FE_Evolution() { }
233 };
234 
235 
236 int main(int argc, char *argv[])
237 {
238  // 1. Initialize MPI.
239  int num_procs, myid;
240  MPI_Init(&argc, &argv);
241  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
242  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
243 
244  // 2. Parse command-line options.
245  problem = 0;
246  optimizer_type = 1;
247  const char *mesh_file = "../../data/periodic-hexagon.mesh";
248  int ser_ref_levels = 2;
249  int par_ref_levels = 0;
250  int order = 3;
251  int ode_solver_type = 3;
252  double t_final = 1.0;
253  double dt = 0.01;
254  bool visualization = true;
255  bool visit = false;
256  bool binary = false;
257  int vis_steps = 5;
258 
259  int precision = 8;
260  cout.precision(precision);
261 
262  OptionsParser args(argc, argv);
263  args.AddOption(&mesh_file, "-m", "--mesh",
264  "Mesh file to use.");
265  args.AddOption(&problem, "-p", "--problem",
266  "Problem setup to use. See options in velocity_function().");
267  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
268  "Number of times to refine the mesh uniformly in serial.");
269  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
270  "Number of times to refine the mesh uniformly in parallel.");
271  args.AddOption(&order, "-o", "--order",
272  "Order (degree) of the finite elements.");
273  args.AddOption(&optimizer_type, "-opt", "--optimizer",
274  "Nonlinear optimizer: 1 - SLBQP,\n\t"
275  " 2 - HIOP.");
276  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
277  "ODE solver: 1 - Forward Euler,\n\t"
278  " 2 - RK2 SSP, 3 - RK3 SSP.");
279  args.AddOption(&t_final, "-tf", "--t-final",
280  "Final time; start time is 0.");
281  args.AddOption(&dt, "-dt", "--time-step",
282  "Time step.");
283  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
284  "--no-visualization",
285  "Enable or disable GLVis visualization.");
286  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
287  "--no-visit-datafiles",
288  "Save data files for VisIt (visit.llnl.gov) visualization.");
289  args.AddOption(&binary, "-binary", "--binary-datafiles", "-ascii",
290  "--ascii-datafiles",
291  "Use binary (Sidre) or ascii format for VisIt data files.");
292  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
293  "Visualize every n-th timestep.");
294  args.Parse();
295  if (!args.Good())
296  {
297  if (myid == 0) { args.PrintUsage(cout); }
298  MPI_Finalize();
299  return 1;
300  }
301  if (myid == 0) { args.PrintOptions(cout); }
302 
303  // 3. Read the serial mesh from the given mesh file on all processors. We can
304  // handle geometrically periodic meshes in this code.
305  Mesh *mesh = new Mesh(mesh_file, 1, 1);
306  int dim = mesh->Dimension();
307 
308  // 4. Define the ODE solver used for time integration. Several explicit
309  // Runge-Kutta methods are available.
310  ODESolver *ode_solver = NULL;
311  switch (ode_solver_type)
312  {
313  case 1: ode_solver = new ForwardEulerSolver; break;
314  case 2: ode_solver = new RK2Solver(1.0); break;
315  case 3: ode_solver = new RK3SSPSolver; break;
316  case 4: ode_solver = new RK4Solver; break;
317  case 6: ode_solver = new RK6Solver; break;
318  default:
319  if (myid == 0)
320  {
321  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
322  }
323  delete mesh;
324  MPI_Finalize();
325  return 3;
326  }
327 
328  // 5. Refine the mesh in serial to increase the resolution. In this example
329  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
330  // a command-line parameter. If the mesh is of NURBS type, we convert it
331  // to a (piecewise-polynomial) high-order mesh.
332  for (int lev = 0; lev < ser_ref_levels; lev++)
333  {
334  mesh->UniformRefinement();
335  }
336  if (mesh->NURBSext)
337  {
338  mesh->SetCurvature(max(order, 1));
339  }
340  mesh->GetBoundingBox(bb_min, bb_max, max(order, 1));
341 
342  // 6. Define the parallel mesh by a partitioning of the serial mesh. Refine
343  // this mesh further in parallel to increase the resolution. Once the
344  // parallel mesh is defined, the serial mesh can be deleted.
345  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
346  delete mesh;
347  for (int lev = 0; lev < par_ref_levels; lev++)
348  {
349  pmesh->UniformRefinement();
350  }
351 
352  // 7. Define the parallel discontinuous DG finite element space on the
353  // parallel refined mesh of the given polynomial order.
354  DG_FECollection fec(order, dim, BasisType::Positive);
355  ParFiniteElementSpace *fes = new ParFiniteElementSpace(pmesh, &fec);
356 
357  HYPRE_Int global_vSize = fes->GlobalTrueVSize();
358  if (myid == 0)
359  {
360  cout << "Number of unknowns: " << global_vSize << endl;
361  }
362 
363  // 8. Set up and assemble the parallel bilinear and linear forms (and the
364  // parallel hypre matrices) corresponding to the DG discretization. The
365  // DGTraceIntegrator involves integrals over mesh interior faces.
369 
370  ParBilinearForm *m = new ParBilinearForm(fes);
372  ParBilinearForm *k = new ParBilinearForm(fes);
373  k->AddDomainIntegrator(new ConvectionIntegrator(velocity, -1.0));
375  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
377  new TransposeIntegrator(new DGTraceIntegrator(velocity, 1.0, -0.5)));
378 
379  ParLinearForm *b = new ParLinearForm(fes);
381  new BoundaryFlowIntegrator(inflow, velocity, -1.0, -0.5));
382 
383  m->Assemble();
384  m->Finalize();
385  int skip_zeros = 0;
386  k->Assemble(skip_zeros);
387  k->Finalize(skip_zeros);
388  b->Assemble();
389 
393 
394  // 9. Define the initial conditions, save the corresponding grid function to
395  // a file and (optionally) save data in the VisIt format and initialize
396  // GLVis visualization.
397  ParGridFunction *u = new ParGridFunction(fes);
398  u->ProjectCoefficient(u0);
399  HypreParVector *U = u->GetTrueDofs();
400 
401  {
402  ostringstream mesh_name, sol_name;
403  mesh_name << "ex9-mesh." << setfill('0') << setw(6) << myid;
404  sol_name << "ex9-init." << setfill('0') << setw(6) << myid;
405  ofstream omesh(mesh_name.str().c_str());
406  omesh.precision(precision);
407  pmesh->Print(omesh);
408  ofstream osol(sol_name.str().c_str());
409  osol.precision(precision);
410  u->Save(osol);
411  }
412 
413  // Create data collection for solution output: either VisItDataCollection for
414  // ascii data files, or SidreDataCollection for binary data files.
415  DataCollection *dc = NULL;
416  if (visit)
417  {
418  if (binary)
419  {
420 #ifdef MFEM_USE_SIDRE
421  dc = new SidreDataCollection("Example9-Parallel", pmesh);
422 #else
423  MFEM_ABORT("Must build with MFEM_USE_SIDRE=YES for binary output.");
424 #endif
425  }
426  else
427  {
428  dc = new VisItDataCollection("Example9-Parallel", pmesh);
429  dc->SetPrecision(precision);
430  // To save the mesh using MFEM's parallel mesh format:
431  // dc->SetFormat(DataCollection::PARALLEL_FORMAT);
432  }
433  dc->RegisterField("solution", u);
434  dc->SetCycle(0);
435  dc->SetTime(0.0);
436  dc->Save();
437  }
438 
439  socketstream sout;
440  if (visualization)
441  {
442  char vishost[] = "localhost";
443  int visport = 19916;
444  sout.open(vishost, visport);
445  if (!sout)
446  {
447  if (myid == 0)
448  cout << "Unable to connect to GLVis server at "
449  << vishost << ':' << visport << endl;
450  visualization = false;
451  if (myid == 0)
452  {
453  cout << "GLVis visualization disabled.\n";
454  }
455  }
456  else
457  {
458  sout << "parallel " << num_procs << " " << myid << "\n";
459  sout.precision(precision);
460  sout << "solution\n" << *pmesh << *u;
461  sout << "pause\n";
462  sout << flush;
463  if (myid == 0)
464  cout << "GLVis visualization paused."
465  << " Press space (in the GLVis window) to resume it.\n";
466  }
467  }
468 
469  Vector M_rowsums(m->Size());
470  m->SpMat().GetRowSums(M_rowsums);
471 
472  // 10. Define the time-dependent evolution operator describing the ODE
473  // right-hand side, and perform time-integration (looping over the time
474  // iterations, ti, with a time-step dt).
475  FE_Evolution adv(*M, *K, *B, *k, M_rowsums);
476 
477  double t = 0.0;
478  adv.SetTime(t);
479  ode_solver->Init(adv);
480 
481  *u = *U;
482 
483  // Compute initial volume.
484  const double vol0_loc = M_rowsums * (*u);
485  double vol0;
486  MPI_Allreduce(&vol0_loc, &vol0, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
487 
488  bool done = false;
489  for (int ti = 0; !done; )
490  {
491  double dt_real = min(dt, t_final - t);
492  adv.SetTimeStep(dt_real);
493  ode_solver->Step(*U, t, dt_real);
494  ti++;
495 
496  done = (t >= t_final - 1e-8*dt);
497 
498  if (done || ti % vis_steps == 0)
499  {
500  if (myid == 0)
501  {
502  cout << "time step: " << ti << ", time: " << t << endl;
503  }
504 
505  // 11. Extract the parallel grid function corresponding to the finite
506  // element approximation U (the local solution on each processor).
507  *u = *U;
508 
509  if (visualization)
510  {
511  sout << "parallel " << num_procs << " " << myid << "\n";
512  sout << "solution\n" << *pmesh << *u << flush;
513  }
514 
515  if (visit)
516  {
517  dc->SetCycle(ti);
518  dc->SetTime(t);
519  dc->Save();
520  }
521  }
522  }
523 
524  // Print the error vs exact solution.
525  const double max_error = u->ComputeMaxError(u0),
526  l1_error = u->ComputeL1Error(u0),
527  l2_error = u->ComputeL2Error(u0);
528  if (myid == 0)
529  {
530  std::cout << "Linf error = " << max_error << endl
531  << "L1 error = " << l1_error << endl
532  << "L2 error = " << l2_error << endl;
533  }
534 
535  // Print error in volume.
536  const double vol_loc = M_rowsums * (*u);
537  double vol;
538  MPI_Allreduce(&vol_loc, &vol, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
539  if (myid == 0)
540  {
541  std::cout << "Vol error = " << vol - vol0 << endl;
542  }
543 
544  // 12. Save the final solution in parallel. This output can be viewed later
545  // using GLVis: "glvis -np <np> -m ex9-mesh -g ex9-final".
546  {
547  *u = *U;
548  ostringstream sol_name;
549  sol_name << "ex9-final." << setfill('0') << setw(6) << myid;
550  ofstream osol(sol_name.str().c_str());
551  osol.precision(precision);
552  u->Save(osol);
553  }
554 
555  // 13. Free the used memory.
556  delete U;
557  delete u;
558  delete B;
559  delete b;
560  delete K;
561  delete k;
562  delete M;
563  delete m;
564  delete fes;
565  delete pmesh;
566  delete ode_solver;
567  delete dc;
568 
569  MPI_Finalize();
570  return 0;
571 }
572 
573 
574 // Implementation of class FE_Evolution
576  const Vector &_b, ParBilinearForm &_pbf,
577  Vector &M_rs)
578  : TimeDependentOperator(_M.Height()),
579  M(_M), K(_K), b(_b), M_solver(M.GetComm()), z(_M.Height()),
580  pbf(_pbf), M_rowsums(M_rs)
581 {
582  M_prec.SetType(HypreSmoother::Jacobi);
583  M_solver.SetPreconditioner(M_prec);
584  M_solver.SetOperator(M);
585 
586  M_solver.iterative_mode = false;
587  M_solver.SetRelTol(1e-9);
588  M_solver.SetAbsTol(0.0);
589  M_solver.SetMaxIter(100);
590  M_solver.SetPrintLevel(0);
591 }
592 
593 void FE_Evolution::Mult(const Vector &x, Vector &y) const
594 {
595  // Get values on the ldofs.
596  ParFiniteElementSpace *pfes = pbf.ParFESpace();
597  ParGridFunction x_gf(pfes);
598  pfes->GetProlongationMatrix()->Mult(x, x_gf);
599 
600  // Compute bounds y_min, y_max for y from from x on the ldofs.
601  const int ldofs = x_gf.Size();
602  Vector y_min(ldofs), y_max(ldofs);
603  x_gf.ExchangeFaceNbrData();
604  Vector &x_nd = x_gf.FaceNbrData();
605  const int *In = pbf.SpMat().GetI(), *Jn = pbf.SpMat().GetJ();
606  for (int i = 0, k = 0; i < ldofs; i++)
607  {
608  double x_i_min = +std::numeric_limits<double>::infinity();
609  double x_i_max = -std::numeric_limits<double>::infinity();
610  for (int end = In[i+1]; k < end; k++)
611  {
612  const int j = Jn[k];
613  const double x_j = (j < ldofs) ? x(j): x_nd(j-ldofs);
614 
615  if (x_j > x_i_max) { x_i_max = x_j; }
616  if (x_j < x_i_min) { x_i_min = x_j; }
617  }
618  y_min(i) = x_i_min;
619  y_max(i) = x_i_max;
620  }
621  for (int i = 0; i < ldofs; i++)
622  {
623  y_min(i) = (y_min(i) - x_gf(i) ) / dt;
624  y_max(i) = (y_max(i) - x_gf(i) ) / dt;
625  }
626  Vector y_min_tdofs(y.Size()), y_max_tdofs(y.Size());
627  // Move the bounds to the tdofs.
628  pfes->GetRestrictionMatrix()->Mult(y_min, y_min_tdofs);
629  pfes->GetRestrictionMatrix()->Mult(y_max, y_max_tdofs);
630 
631  // Compute the high-order solution y = M^{-1} (K x + b) on the tdofs.
632  K.Mult(x, z);
633  z += b;
634  M_solver.Mult(z, y);
635 
636  // The solution y is an increment; it should not introduce new mass.
637  const double mass_y = 0.0;
638 
639  // Perform optimization on the tdofs.
640  Vector y_out(y.Size());
641  const int max_iter = 500;
642  const double rtol = 1.e-7;
643  double atol = 1.e-7;
644 
645  OptimizationSolver* optsolver = NULL;
646  if (optimizer_type == 2)
647  {
648 #ifdef MFEM_USE_HIOP
649  HiopNlpOptimizer *tmp_opt_ptr = new HiopNlpOptimizer(MPI_COMM_WORLD);
650  optsolver = tmp_opt_ptr;
651 #else
652  MFEM_ABORT("MFEM is not built with HiOp support!");
653 #endif
654  }
655  else
656  {
657  SLBQPOptimizer *slbqp = new SLBQPOptimizer(MPI_COMM_WORLD);
658  slbqp->SetBounds(y_min_tdofs, y_max_tdofs);
659  slbqp->SetLinearConstraint(M_rowsums, mass_y);
660  atol = 1.e-15;
661  optsolver = slbqp;
662  }
663 
664  OptimizedTransportProblem ot_prob(*pfes, y, M_rowsums, mass_y,
665  y_min_tdofs, y_max_tdofs);
666  optsolver->SetOptimizationProblem(ot_prob);
667 
668  optsolver->SetMaxIter(max_iter);
669  optsolver->SetAbsTol(atol);
670  optsolver->SetRelTol(rtol);
671  optsolver->SetPrintLevel(0);
672  optsolver->Mult(y, y_out);
673 
674  y = y_out;
675 
676  delete optsolver;
677 }
678 
679 
680 // Velocity coefficient
681 void velocity_function(const Vector &x, Vector &v)
682 {
683  int dim = x.Size();
684 
685  // map to the reference [-1,1] domain
686  Vector X(dim);
687  for (int i = 0; i < dim; i++)
688  {
689  double center = (bb_min[i] + bb_max[i]) * 0.5;
690  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
691  }
692 
693  switch (problem)
694  {
695  case 0:
696  {
697  // Translations in 1D, 2D, and 3D
698  switch (dim)
699  {
700  case 1: v(0) = (invert_velocity) ? -1.0 : 1.0; break;
701  case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.); break;
702  case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
703  break;
704  }
705  break;
706  }
707  case 1:
708  case 2:
709  {
710  // Clockwise rotation in 2D around the origin
711  const double w = M_PI/2;
712  switch (dim)
713  {
714  case 1: v(0) = 1.0; break;
715  case 2: v(0) = w*X(1); v(1) = -w*X(0); break;
716  case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0; break;
717  }
718  break;
719  }
720  case 3:
721  {
722  // Clockwise twisting rotation in 2D around the origin
723  const double w = M_PI/2;
724  double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
725  d = d*d;
726  switch (dim)
727  {
728  case 1: v(0) = 1.0; break;
729  case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0); break;
730  case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0; break;
731  }
732  break;
733  }
734  }
735 }
736 
737 // Initial condition
738 double u0_function(const Vector &x)
739 {
740  int dim = x.Size();
741 
742  // map to the reference [-1,1] domain
743  Vector X(dim);
744  for (int i = 0; i < dim; i++)
745  {
746  double center = (bb_min[i] + bb_max[i]) * 0.5;
747  X(i) = 2 * (x(i) - center) / (bb_max[i] - bb_min[i]);
748  }
749 
750  switch (problem)
751  {
752  case 0:
753  case 1:
754  {
755  switch (dim)
756  {
757  case 1:
758  return (X(0) > -0.15 && X(0) < 0.15) ? 1.0 : 0.0;
759  //return exp(-40.*pow(X(0)-0.0,2));
760  case 2:
761  case 3:
762  {
763  double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
764  if (dim == 3)
765  {
766  const double s = (1. + 0.25*cos(2*M_PI*X(2)));
767  rx *= s;
768  ry *= s;
769  }
770  return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
771  erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
772  }
773  }
774  }
775  case 2:
776  {
777  double x_ = X(0), y_ = X(1), rho, phi;
778  rho = hypot(x_, y_);
779  phi = atan2(y_, x_);
780  return pow(sin(M_PI*rho),2)*sin(3*phi);
781  }
782  case 3:
783  {
784  const double f = M_PI;
785  return sin(f*X(0))*sin(f*X(1));
786  }
787  }
788  return 0.0;
789 }
790 
791 // Inflow boundary condition (zero for the problems considered in this example)
792 double inflow_function(const Vector &x)
793 {
794  switch (problem)
795  {
796  case 0:
797  case 1:
798  case 2:
799  case 3: return 0.0;
800  }
801  return 0.0;
802 }
Vector bb_max
Definition: ex9.cpp:60
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
Conjugate gradient method.
Definition: solvers.hpp:152
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:358
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:890
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
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:485
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:301
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
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
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
Class for parallel linear form.
Definition: plinearform.hpp:26
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
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
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
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4102
Parallel smoothers in hypre.
Definition: hypre.hpp:581
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: ex9p.cpp:228
Vector bb_min
Definition: ex9.cpp:60
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
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
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:548
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
string space
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:132
virtual ~FE_Evolution()
Definition: ex9p.cpp:232
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 const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:339
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
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
Class for parallel bilinear form.
void SetK(HypreParMatrix &_K)
Definition: ex9p.cpp:229
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
Vector data type.
Definition: vector.hpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The classical forward Euler method.
Definition: ode.hpp:99
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
Definition: solvers.hpp:408
Class for parallel meshes.
Definition: pmesh.hpp:32
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:816
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
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)