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