77#ifndef MFEM_USE_SUNDIALS
78#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
100 void QuadratureIntegration(
const Vector &x,
Vector &y)
const override;
104 Vector &yBdot)
const override;
107 void QuadratureSensitivityMult(
const Vector &y,
const Vector &yB,
108 Vector &qbdot)
const override;
114 int *jcurB,
real_t gammaB)
override;
119 ~RobertsTDAOperator()
121 delete adjointMatrix;
132int main(
int argc,
char *argv[])
139 const real_t reltol = 1e-4;
142 cout.precision(precision);
145 args.
AddOption(&dt,
"-dt",
"--time-step",
"Time step.");
172 RobertsTDAOperator adv(3,
p);
181 abstol_v[0] = 1.0e-8;
182 abstol_v[1] = 1.0e-14;
183 abstol_v[2] = 1.0e-6;
199 for (
int i = 0; i < y.
Size(); i++)
201 real_t ww = reltol * abs(y[i]) + abstol_v[i];
202 if (ww <= 0.) {
return -1; }
226 real_t dt_real = max(dt, t_final - t);
227 cvodes->
Step(
u, t, dt_real);
229 done = (t >= t_final - 1e-8*dt);
237 cout <<
"Final Solution: " << t << endl;
241 cout <<
" Final Quadrature " << endl;
261 real_t dt_real = max(dt, t - TBout1);
262 cvodes->
StepB(w, t, dt_real);
263 cout <<
"t: " << t << endl;
264 cout <<
"w:" << endl;
268 cout <<
"u:" << endl;
272 dt_real = max(dt, t - 0.);
273 cvodes->
StepB(w, t, dt_real);
274 cout <<
"t: " << t << endl;
275 cout <<
"w:" << endl;
279 cout <<
"u:" << endl;
284 cout <<
"dG/dp:" << endl;
296void RobertsTDAOperator::Mult(
const Vector &x,
Vector &y)
const
298 y[0] = -p_[0]*x[0] + p_[1]*x[1]*x[2];
299 y[2] = p_[2]*x[1]*x[1];
304void RobertsTDAOperator::QuadratureIntegration(
const Vector &y,
311void RobertsTDAOperator::AdjointRateMult(
const Vector &y,
Vector & yB,
314 real_t l21 = (yB[1]-yB[0]);
315 real_t l32 = (yB[2]-yB[1]);
319 yBdot[0] = -p1 * l21;
320 yBdot[1] = p2 * y[2] * l21 - 2. * p3 * y[1] * l32;
321 yBdot[2] = p2 * y[1] * l21 - 1.0;
325void RobertsTDAOperator::QuadratureSensitivityMult(
const Vector &y,
329 real_t l21 = (yB[1]-yB[0]);
330 real_t l32 = (yB[2]-yB[1]);
333 qBdot[0] = y[0] * l21;
334 qBdot[1] = -y23 * l21;
335 qBdot[2] = y[1]*y[1]*l32;
339int RobertsTDAOperator::SUNImplicitSetupB(
const real_t t,
const Vector &y,
341 const Vector &fyB,
int jokB,
342 int *jcurB,
real_t gammaB)
346 delete adjointMatrix;
348 for (
int j = 0; j < y.
Size(); j++)
354 AdjointRateMult(y, yBone, JacBj);
356 for (
int i = 0; i < y.
Size(); i++)
358 adjointMatrix->
Set(i,j, (i == j ? 1.0 : 0.) - gammaB * JacBj[i]);
370int RobertsTDAOperator::SUNImplicitSolveB(
Vector &x,
const Vector &
b,
375 adjointSolver.
Mult(
b, x);
void EvalQuadIntegrationB(double t, Vector &dG_dp)
Evaluate Quadrature solution.
void EvalQuadIntegration(double t, Vector &q)
Evaluate Quadrature.
void SetMaxNStepsB(int mxstepsB)
Set the maximum number of backward steps.
void Step(Vector &x, double &t, double &dt) override
void InitB(TimeDependentAdjointOperator &f_)
Initialize the adjoint problem.
void GetForwardSolution(double tB, mfem::Vector &yy)
Get Interpolated Forward solution y at backward integration time tB.
void SetWFTolerances(EWTFunction func)
Set multiplicative error weights.
void Init(TimeDependentAdjointOperator &f_)
virtual void StepB(Vector &w, double &t, double &dt)
Solve one adjoint time step.
void InitQuadIntegrationB(mfem::Vector &qB0, double reltolQB=1e-3, double abstolQB=1e-8)
Initialize Quadrature Integration (Adjoint)
void InitAdjointSolve(int steps, int interpolation)
Initialize Adjoint.
void InitQuadIntegration(mfem::Vector &q0, double reltolQ=1e-3, double abstolQ=1e-8)
Interface to the CVODE library – linear multi-step methods.
void PrintInfo() const
Print various CVODE statistics.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
void SetSVtolerances(double reltol, Vector abstol)
Set the scalar relative and vector of absolute tolerances.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the GMRES method.
void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
void SetRelTol(real_t rtol)
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.
void Finalize(int skip_zeros=1) override
Finalize the matrix initialization, switching the storage format from LIL to CSR.
void Set(const int i, const int j, const real_t val)
virtual void SetTime(const real_t t_)
Set the current time.
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
int Size() const
Returns the size of the vector.
real_t u(const Vector &xvec)
real_t p(const Vector &x, real_t t)