73 #ifndef MFEM_USE_SUNDIALS 74 #error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES 84 AdvDiffSUNDIALS(
int ydot_dim,
int ybdot_dim,
Vector p,
88 ess_tdof_list(ess_tdof),
91 M_solver(fes->GetComm())
95 cout <<
"Essential tdofs: " << endl;
96 ess_tdof_list.Print();
100 m->Assemble(skip_zeros);
101 m->Finalize(skip_zeros);
113 k->Assemble(skip_zeros);
114 k->Finalize(skip_zeros);
119 k1->Assemble(skip_zeros);
120 k1->Finalize(skip_zeros);
122 M = m->ParallelAssemble();
126 K = k->ParallelAssemble();
130 K_adj = k1->ParallelAssemble();
134 M_prec.SetType(HypreSmoother::Jacobi);
135 M_solver.SetPreconditioner(M_prec);
136 M_solver.SetOperator(*M);
138 M_solver.SetRelTol(1e-14);
139 M_solver.SetAbsTol(0.0);
140 M_solver.SetMaxIter(1000);
141 M_solver.SetPrintLevel(0);
147 virtual void AdjointRateMult(
const Vector &y,
Vector &yB,
150 virtual int SUNImplicitSetup(
const Vector &y,
151 const Vector &fy,
int jok,
int *jcur,
154 virtual int SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol);
156 virtual void QuadratureSensitivityMult(
const Vector &y,
const Vector &yB,
199 return x[0]*(2. - x[0])*exp(2.*x[0]);
202 int main(
int argc,
char *argv[])
205 Mpi::Init(argc, argv);
206 int myid = Mpi::WorldRank();
210 int ser_ref_levels = 0;
211 int par_ref_levels = 0;
212 double t_final = 2.5;
215 bool step_mode =
true;
218 cout.precision(precision);
222 args.
AddOption(&mx,
"-m",
"--mx",
"The number of mesh elements in the x-dir");
223 args.
AddOption(&ser_ref_levels,
"-r",
"--refine",
224 "Number of times to refine the mesh uniformly.");
225 args.
AddOption(&step_mode,
"-a",
"--adams",
"-no-a",
"--no-adams",
226 "A switch to toggle between CV_ADAMS, and CV_BDF stepping modes");
227 args.
AddOption(&t_final,
"-tf",
"--t-final",
228 "Final time; start time is 0.");
229 args.
AddOption(&dt,
"-dt",
"--time-step",
248 Mesh mesh = Mesh::MakeCartesian1D(mx+1, 2.);
254 for (
int lev = 0; lev < ser_ref_levels; lev++)
261 for (
int lev = 0; lev < par_ref_levels; lev++)
273 cout <<
"Number of unknowns: " << global_vSize << endl;
286 u.ProjectCoefficient(u0);
288 cout <<
"Init u: " << endl;
297 essential_attr[0] = 1;
298 essential_attr[1] = 1;
299 fes->GetEssentialTrueDofs(essential_attr, ess_tdof_list);
302 AdvDiffSUNDIALS adv(U->
Size(), U->
Size(),
p, fes, ess_tdof_list);
310 step_mode ? CV_ADAMS : CV_BDF);
316 double reltol = 1e-8, abstol = 1e-6;
320 int checkpoint_steps = 50;
326 for (
int ti = 0; !done; )
328 double dt_real = max(dt, t_final -
t);
329 cvodes->
Step(*U,
t, dt_real);
332 done = (
t >= t_final - 1e-8*dt);
334 if (done && myid == 0)
343 cout <<
"Final Solution: " <<
t << endl;
346 cout <<
"u (" << myid <<
"):" << endl;
349 MPI_Barrier(MPI_COMM_WORLD);
361 cout <<
"g: " << g << endl;
381 double dt_real = max(dt,
t);
382 cvodes->
StepB(*V,
t, dt_real);
385 cout <<
"t: " <<
t << endl;
388 cout <<
"v (" << myid <<
"):" << endl;
391 MPI_Barrier(MPI_COMM_WORLD);
396 MPI_Barrier(MPI_COMM_WORLD);
399 cout <<
"sensitivity:" << endl;
420 x1.SetSubVector(ess_tdof_list, 0.0);
429 int AdvDiffSUNDIALS::SUNImplicitSetup(
const Vector &y,
431 int jok,
int *jcur,
double gamma)
438 Mf =
Add(1., *M, -gamma, *K);
445 int AdvDiffSUNDIALS::SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol)
452 prec.
SetType(HypreSmoother::Jacobi);
453 solver.SetPreconditioner(prec);
454 solver.SetOperator(*Mf);
455 solver.SetRelTol(1E-14);
456 solver.SetMaxIter(1000);
464 void AdvDiffSUNDIALS::AdjointRateMult(
const Vector &y,
Vector & yB,
472 M_solver.Mult(z, yBdot);
476 void AdvDiffSUNDIALS::QuadratureSensitivityMult(
const Vector &y,
503 Vector p2vec(pfes->GetParMesh()->SpaceDimension()); p2vec = 1.;
517 double dp1_result =
InnerProduct(pfes->GetComm(), yB, b1);
518 double dp2_result =
InnerProduct(pfes->GetComm(), yB, b2);
520 qBdot[0] = -dp1_result;
521 qBdot[1] = -dp2_result;
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Class for domain integration L(v) := (f, v)
Conjugate gradient method.
virtual void StepB(Vector &w, double &t, double &dt)
Solve one adjoint time step.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
void PrintUsage(std::ostream &out) const
Print the usage message.
Vector coefficient that is constant in space and time.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Size() const
Returns the size of the vector.
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
virtual void Step(Vector &x, double &t, double &dt)
void UseSundialsLinearSolverB()
Use built in SUNDIALS Newton solver.
void SetSStolerancesB(double reltol, double abstol)
Tolerance specification functions for the adjoint problem.
void InitAdjointSolve(int steps, int interpolation)
Initialize Adjoint.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
int main(int argc, char *argv[])
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void InitQuadIntegrationB(mfem::Vector &qB0, double reltolQB=1e-3, double abstolQB=1e-8)
Initialize Quadrature Integration (Adjoint)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
ParMesh * GetParMesh() const
Parallel smoothers in hypre.
Wrapper for hypre's parallel vector class.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
double p(const Vector &x, double t)
void PrintInfo() const
Print various CVODE statistics.
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...
int SpaceDimension() const
Dimension of the physical space containing the mesh.
double InnerProduct(HypreParVector *x, HypreParVector *y)
void InitB(TimeDependentAdjointOperator &f_)
Initialize the adjoint problem.
int Size() const
Return the logical size of the array.
void Clear()
Clear the contents of the Mesh.
A general function coefficient.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Arbitrary order H1-conforming (continuous) finite elements.
void Init(TimeDependentAdjointOperator &f_)
double u(const Vector &xvec)
Class for parallel grid function.
Wrapper for hypre's ParCSR matrix class.
double u_init(const Vector &x)
Class for parallel meshes.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
void EvalQuadIntegrationB(double t, Vector &dG_dp)
Evaluate Quadrature solution.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
void SetMaxNSteps(int steps)
Set the maximum number of time steps.