77 #ifndef MFEM_USE_SUNDIALS 78 #error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES 100 virtual void QuadratureIntegration(
const Vector &x,
Vector &y)
const;
103 virtual void AdjointRateMult(
const Vector &y,
Vector &yB,
107 virtual void QuadratureSensitivityMult(
const Vector &y,
const Vector &yB,
112 virtual int SUNImplicitSetupB(
const double t,
const Vector &y,
114 int *jcurB,
double gammaB);
117 virtual int SUNImplicitSolveB(
Vector &x,
const Vector &
b,
double tol);
119 ~RobertsTDAOperator()
121 delete adjointMatrix;
132 int main(
int argc,
char *argv[])
135 double t_final = 4e7;
139 const double 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 double ww = reltol * abs(y[i]) + abstol_v[i];
202 if (ww <= 0.) {
return -1; }
224 for (
int ti = 0; !done; )
226 double dt_real = max(dt, t_final -
t);
227 cvodes->
Step(
u,
t, dt_real);
230 done = (
t >= t_final - 1e-8*dt);
238 cout <<
"Final Solution: " <<
t << endl;
242 cout <<
" Final Quadrature " << endl;
262 double dt_real = max(dt,
t - TBout1);
263 cvodes->
StepB(w,
t, dt_real);
264 cout <<
"t: " <<
t << endl;
265 cout <<
"w:" << endl;
269 cout <<
"u:" << endl;
273 dt_real = max(dt,
t - 0.);
274 cvodes->
StepB(w,
t, dt_real);
275 cout <<
"t: " <<
t << endl;
276 cout <<
"w:" << endl;
280 cout <<
"u:" << endl;
285 cout <<
"dG/dp:" << endl;
299 y[0] = -p_[0]*x[0] + p_[1]*x[1]*x[2];
300 y[2] = p_[2]*x[1]*x[1];
305 void RobertsTDAOperator::QuadratureIntegration(
const Vector &y,
312 void RobertsTDAOperator::AdjointRateMult(
const Vector &y,
Vector & yB,
315 double l21 = (yB[1]-yB[0]);
316 double l32 = (yB[2]-yB[1]);
320 yBdot[0] = -p1 * l21;
321 yBdot[1] = p2 * y[2] * l21 - 2. * p3 * y[1] * l32;
322 yBdot[2] = p2 * y[1] * l21 - 1.0;
326 void RobertsTDAOperator::QuadratureSensitivityMult(
const Vector &y,
330 double l21 = (yB[1]-yB[0]);
331 double l32 = (yB[2]-yB[1]);
332 double y23 = y[1] * y[2];
334 qBdot[0] = y[0] * l21;
335 qBdot[1] = -y23 * l21;
336 qBdot[2] = y[1]*y[1]*l32;
340 int RobertsTDAOperator::SUNImplicitSetupB(
const double t,
const Vector &y,
342 const Vector &fyB,
int jokB,
343 int *jcurB,
double gammaB)
347 delete adjointMatrix;
349 for (
int j = 0; j < y.
Size(); j++)
355 AdjointRateMult(y, yBone, JacBj);
357 for (
int i = 0; i < y.
Size(); i++)
359 adjointMatrix->Set(i,j, (i == j ? 1.0 : 0.) - gammaB * JacBj[i]);
364 adjointMatrix->Finalize();
365 adjointSolver.SetOperator(*adjointMatrix);
371 int RobertsTDAOperator::SUNImplicitSolveB(
Vector &x,
const Vector &
b,
375 adjointSolver.SetRelTol(1e-14);
376 adjointSolver.Mult(
b, x);
void GetForwardSolution(double tB, mfem::Vector &yy)
Get Interpolated Forward solution y at backward integration time tB.
virtual void StepB(Vector &w, double &t, double &dt)
Solve one adjoint time step.
void PrintOptions(std::ostream &out) const
Print the options.
void PrintUsage(std::ostream &out) const
Print the usage message.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
void SetSVtolerances(double reltol, Vector abstol)
Set the scalar relative and vector of absolute tolerances.
int Size() const
Returns the size of the vector.
void EvalQuadIntegration(double t, Vector &q)
Evaluate Quadrature.
bool Good() const
Return true if the command line options were parsed successfully.
virtual void Step(Vector &x, double &t, double &dt)
void InitAdjointSolve(int steps, int interpolation)
Initialize Adjoint.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Interface to the CVODE library – linear multi-step methods.
void InitQuadIntegration(mfem::Vector &q0, double reltolQ=1e-3, double abstolQ=1e-8)
void InitQuadIntegrationB(mfem::Vector &qB0, double reltolQB=1e-3, double abstolQB=1e-8)
Initialize Quadrature Integration (Adjoint)
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...
void InitB(TimeDependentAdjointOperator &f_)
Initialize the adjoint problem.
void SetWFTolerances(EWTFunction func)
Set multiplicative error weights.
int main(int argc, char *argv[])
void Init(TimeDependentAdjointOperator &f_)
double u(const Vector &xvec)
void SetMaxNStepsB(int mxstepsB)
Set the maximum number of backward steps.
void EvalQuadIntegrationB(double t, Vector &dG_dp)
Evaluate Quadrature solution.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.