56#error This example requires that MFEM is built with MFEM_USE_HIOP=YES
83class LinearScaleOperator :
public Operator
99 for (
int i = 0; i <
width; i++) { grad(0, i) = w_glob(i); }
106 const double loc_res = w * x_loc;
107 MPI_Allreduce(&loc_res, &y(0), 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
118class TanhSumOperator :
public Operator
130 double sum_loc = x.
Sum();
131 MPI_Allreduce(&sum_loc, &y(0), 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
132 y(0) = std::tanh(y(0));
137 double sum_loc = x.
Sum();
139 MPI_Allreduce(&sum_loc, &dtanh, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
140 dtanh = 1.0 - pow(std::tanh(dtanh), 2);
142 for (
int i = 0; i <
width; i++) { grad(0, i) = dtanh; }
158 Vector massvec, d_lo, d_hi;
159 const LinearScaleOperator LSoper;
160 const TanhSumOperator TSoper;
167 x_HO(xho), massvec(1), d_lo(1), d_hi(1),
175 double lsums[2], gsums[2];
176 lsums[0] = xmin.
Sum();
177 lsums[1] = xmax.
Sum();
178 MPI_Allreduce(lsums, gsums, 2, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
179 d_lo(0) = std::tanh(gsums[0]);
180 d_hi(0) = std::tanh(gsums[1]);
181 MFEM_ASSERT(d_lo(0) < d_hi(0),
182 "The bounds produce an infeasible optimization problem");
188 virtual double CalcObjective(
const Vector &x)
const
190 double loc_res = 0.0;
193 const double d = x(i) - x_HO(i);
198 MPI_Allreduce(&loc_res, &res, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
202 virtual void CalcObjectiveGrad(
const Vector &x,
Vector &grad)
const
204 for (
int i = 0; i <
input_size; i++) { grad(i) = x(i) - x_HO(i); }
232 void SetTimeStep(
double dt_) { dt = dt_; }
236 virtual ~FE_Evolution() { }
240int main(
int argc,
char *argv[])
251 const char *mesh_file =
"../../data/periodic-hexagon.mesh";
252 int ser_ref_levels = 2;
253 int par_ref_levels = 0;
255 int ode_solver_type = 3;
256 double t_final = 1.0;
258 bool visualization =
true;
264 cout.precision(precision);
267 args.
AddOption(&mesh_file,
"-m",
"--mesh",
268 "Mesh file to use.");
270 "Problem setup to use. See options in velocity_function().");
271 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
272 "Number of times to refine the mesh uniformly in serial.");
273 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
274 "Number of times to refine the mesh uniformly in parallel.");
276 "Order (degree) of the finite elements.");
278 "Nonlinear optimizer: 1 - SLBQP,\n\t"
280 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
281 "ODE solver: 1 - Forward Euler,\n\t"
282 " 2 - RK2 SSP, 3 - RK3 SSP.");
283 args.
AddOption(&t_final,
"-tf",
"--t-final",
284 "Final time; start time is 0.");
285 args.
AddOption(&dt,
"-dt",
"--time-step",
287 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
288 "--no-visualization",
289 "Enable or disable GLVis visualization.");
290 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
291 "--no-visit-datafiles",
292 "Save data files for VisIt (visit.llnl.gov) visualization.");
293 args.
AddOption(&binary,
"-binary",
"--binary-datafiles",
"-ascii",
295 "Use binary (Sidre) or ascii format for VisIt data files.");
296 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
297 "Visualize every n-th timestep.");
308 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
314 switch (ode_solver_type)
317 case 2: ode_solver =
new RK2Solver(1.0);
break;
319 case 4: ode_solver =
new RK4Solver;
break;
320 case 6: ode_solver =
new RK6Solver;
break;
324 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
334 for (
int lev = 0; lev < ser_ref_levels; lev++)
349 for (
int lev = 0; lev < par_ref_levels; lev++)
362 cout <<
"Number of unknowns: " << global_vSize << endl;
382 b->AddBdrFaceIntegrator(
400 u->ProjectCoefficient(u0);
404 ostringstream mesh_name, sol_name;
405 mesh_name <<
"ex9-mesh." << setfill(
'0') << setw(6) << myid;
406 sol_name <<
"ex9-init." << setfill(
'0') << setw(6) << myid;
407 ofstream omesh(mesh_name.str().c_str());
408 omesh.precision(precision);
410 ofstream osol(sol_name.str().c_str());
411 osol.precision(precision);
425 MFEM_ABORT(
"Must build with MFEM_USE_SIDRE=YES for binary output.");
450 cout <<
"Unable to connect to GLVis server at "
452 visualization =
false;
455 cout <<
"GLVis visualization disabled.\n";
460 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
461 sout.precision(precision);
462 sout <<
"solution\n" << *pmesh << *
u;
466 cout <<
"GLVis visualization paused."
467 <<
" Press space (in the GLVis window) to resume it.\n";
477 FE_Evolution adv(*M, *K, *B, *k, M_rowsums);
481 ode_solver->
Init(adv);
486 const double vol0_loc = M_rowsums * (*u);
488 MPI_Allreduce(&vol0_loc, &vol0, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
491 for (
int ti = 0; !done; )
493 double dt_real = min(dt, t_final -
t);
494 adv.SetTimeStep(dt_real);
495 ode_solver->
Step(*U,
t, dt_real);
498 done = (
t >= t_final - 1e-8*dt);
500 if (done || ti % vis_steps == 0)
504 cout <<
"time step: " << ti <<
", time: " <<
t << endl;
513 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
514 sout <<
"solution\n" << *pmesh << *
u << flush;
527 const double max_error =
u->ComputeMaxError(u0),
528 l1_error =
u->ComputeL1Error(u0),
529 l2_error =
u->ComputeL2Error(u0);
532 std::cout <<
"Linf error = " << max_error << endl
533 <<
"L1 error = " << l1_error << endl
534 <<
"L2 error = " << l2_error << endl;
538 const double vol_loc = M_rowsums * (*u);
540 MPI_Allreduce(&vol_loc, &vol, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
543 std::cout <<
"Vol error = " << vol - vol0 << endl;
550 ostringstream sol_name;
551 sol_name <<
"ex9-final." << setfill(
'0') << setw(6) << myid;
552 ofstream osol(sol_name.str().c_str());
553 osol.precision(precision);
580 M(M_), K(K_),
b(b_), M_solver(M.GetComm()), z(M_.Height()),
581 pbf(pbf_), M_rowsums(M_rs)
584 M_solver.SetPreconditioner(M_prec);
585 M_solver.SetOperator(M);
587 M_solver.iterative_mode =
false;
588 M_solver.SetRelTol(1e-9);
589 M_solver.SetAbsTol(0.0);
590 M_solver.SetMaxIter(100);
591 M_solver.SetPrintLevel(0);
594void FE_Evolution::Mult(
const Vector &x,
Vector &y)
const
602 const int ldofs = x_gf.Size();
603 Vector y_min(ldofs), y_max(ldofs);
604 x_gf.ExchangeFaceNbrData();
605 Vector &x_nd = x_gf.FaceNbrData();
607 for (
int i = 0, k = 0; i < ldofs; i++)
609 double x_i_min = +std::numeric_limits<double>::infinity();
610 double x_i_max = -std::numeric_limits<double>::infinity();
611 for (
int end = In[i+1]; k < end; k++)
614 const double x_j = (j < ldofs) ? x(j): x_nd(j-ldofs);
616 if (x_j > x_i_max) { x_i_max = x_j; }
617 if (x_j < x_i_min) { x_i_min = x_j; }
622 for (
int i = 0; i < ldofs; i++)
624 y_min(i) = (y_min(i) - x_gf(i) ) / dt;
625 y_max(i) = (y_max(i) - x_gf(i) ) / dt;
638 const double mass_y = 0.0;
642 const int max_iter = 500;
643 const double rtol = 1.e-7;
651 optsolver = tmp_opt_ptr;
653 MFEM_ABORT(
"MFEM is not built with HiOp support!");
659 slbqp->
SetBounds(y_min_tdofs, y_max_tdofs);
665 OptimizedTransportProblem ot_prob(*pfes, y, M_rowsums, mass_y,
666 y_min_tdofs, y_max_tdofs);
673 optsolver->
Mult(y, y_out);
688 for (
int i = 0; i <
dim; i++)
702 case 2: v(0) = sqrt(2./3.); v(1) = sqrt(1./3.);
break;
703 case 3: v(0) = sqrt(3./6.); v(1) = sqrt(2./6.); v(2) = sqrt(1./6.);
712 const double w = M_PI/2;
715 case 1: v(0) = 1.0;
break;
716 case 2: v(0) = w*X(1); v(1) = -w*X(0);
break;
717 case 3: v(0) = w*X(1); v(1) = -w*X(0); v(2) = 0.0;
break;
724 const double w = M_PI/2;
725 double d = max((X(0)+1.)*(1.-X(0)),0.) * max((X(1)+1.)*(1.-X(1)),0.);
729 case 1: v(0) = 1.0;
break;
730 case 2: v(0) = d*w*X(1); v(1) = -d*w*X(0);
break;
731 case 3: v(0) = d*w*X(1); v(1) = -d*w*X(0); v(2) = 0.0;
break;
745 for (
int i = 0; i <
dim; i++)
759 return (X(0) > -0.15 && X(0) < 0.15) ? 1.0 : 0.0;
764 double rx = 0.45, ry = 0.25, cx = 0., cy = -0.2, w = 10.;
767 const double s = (1. + 0.25*cos(2*M_PI*X(2)));
771 return ( erfc(w*(X(0)-cx-rx))*erfc(-w*(X(0)-cx+rx)) *
772 erfc(w*(X(1)-cy-ry))*erfc(-w*(X(1)-cy+ry)) )/16;
778 double x_ = X(0), y_ = X(1), rho, phi;
781 return pow(sin(M_PI*rho),2)*sin(3*phi);
785 const double f = M_PI;
786 return sin(
f*X(0))*sin(
f*X(1));
@ Positive
Bernstein polynomials.
Conjugate gradient method.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
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.
The classical forward Euler method.
A general function coefficient.
Adapts the HiOp functionality to the MFEM OptimizationSolver interface.
Wrapper for hypre's ParCSR matrix class.
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.
Wrapper for hypre's parallel vector class.
Parallel smoothers in hypre.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
void SetRelTol(real_t rtol)
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
void SetAbsTol(real_t atol)
Arbitrary order "L2-conforming" discontinuous finite elements.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
int Dimension() const
Dimension of the reference space used within the elements.
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.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
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)
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
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].
int width
Dimension of the input / number of columns in the matrix.
Operator(int s=0)
Construct a square Operator with given size s (default 0).
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).
void SetSolutionBounds(const Vector &xl, const Vector &xh)
void SetEqualityConstraint(const Vector &c)
void SetInequalityConstraint(const Vector &dl, const Vector &dh)
Abstract solver for OptimizationProblems.
virtual void Mult(const Vector &xt, Vector &x) const =0
Operator application: y=A(x).
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
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...
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
HYPRE_BigInt GlobalTrueVSize() const
const Operator * GetProlongationMatrix() const override
The returned Operator is owned by the FiniteElementSpace.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
Class for parallel grid function.
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Third-order, strong stability preserving (SSP) Runge-Kutta method.
The classical explicit forth-order Runge-Kutta method, RK4.
void SetBounds(const Vector &lo_, const Vector &hi_)
void SetLinearConstraint(const Vector &w_, real_t a_)
Data collection with Sidre routines following the Conduit mesh blueprint specification.
void GetRowSums(Vector &x) const
For all i compute .
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
int * GetJ()
Return the array J.
int * GetI()
Return the array I.
Base abstract class for first order time dependent operators.
virtual void SetTime(const real_t t_)
Set the current time.
A general vector function coefficient.
int Size() const
Returns the size of the vector.
real_t Sum() const
Return the sum of the vector entries.
Data collection with VisIt I/O routines.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
void velocity_function(const Vector &x, Vector &v)
double u0_function(const Vector &x)
double inflow_function(const Vector &x)
real_t u(const Vector &xvec)
std::function< real_t(const Vector &)> f(real_t mass_coeff)