Quandary
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | List of all members
ImplMidpoint Class Reference

Implicit midpoint rule time integration scheme. More...

#include <timestepper.hpp>

Inheritance diagram for ImplMidpoint:
Inheritance graph
[legend]
Collaboration diagram for ImplMidpoint:
Collaboration graph
[legend]

Public Member Functions

 ImplMidpoint (MasterEq *mastereq_, int ntime_, double total_time_, LinearSolverType linsolve_type_, int linsolve_maxiter_, Output *output_, bool storeFWD_)
 Constructor for implicit midpoint scheme.
 
 ~ImplMidpoint ()
 
virtual void evolveFWD (const double tstart, const double tstop, Vec x)
 Evolves state forward using implicit midpoint rule.
 
virtual void evolveBWD (const double tstart, const double tstop, const Vec x_stop, Vec x_adj, Vec grad, bool compute_gradient)
 Evolves adjoint backward using implicit midpoint rule and adds to reduced gradient.
 
int NeumannSolve (Mat A, Vec b, Vec x, double alpha, bool transpose)
 Solves (I - alpha*A) * x = b using Neumann iterations.
 
- Public Member Functions inherited from TimeStepper
 TimeStepper ()
 
 TimeStepper (MasterEq *mastereq_, int ntime_, double total_time_, Output *output_, bool storeFWD_)
 Constructor for time stepper.
 
virtual ~TimeStepper ()
 
Vec getState (size_t tindex)
 Retrieves stored state at a specific time index.
 
void setEvalLeakage (bool flag)
 
void setEvalWeightedCost (bool flag, double width)
 
void setEvalDPDM (bool flag)
 
void setEvalEnergy (bool flag)
 
double getLeakageIntegral ()
 
double getWeightedCostIntegral ()
 
double getEnergyIntegral ()
 
double getDPDMIntegral ()
 
Vec solveODE (int initid, Vec rho_t0)
 Solves the ODE forward in time.
 
void solveAdjointODE (Vec rho_t0_bar, Vec finalstate, double Jbar_leakage, double Jbar_weightedcost, double Jbar_dpdm, double Jbar_energy)
 Solves the adjoint ODE backward in time.
 
double evalLeakage (const Vec x)
 Evaluates leakage into guard levels.
 
void evalLeakage_diff (const Vec x, Vec xbar, double Jbar)
 Computes derivative of leakage term.
 
double evalWeightedCost (double time, const Vec x)
 Evaluates the weighted cost function term.
 
void evalWeightedCost_diff (double time, const Vec x, Vec xbar, double Jbar)
 Computes derivative of weighted cost function term.
 
double evalDpDm (Vec x, Vec xm1, Vec xm2)
 Evaluates second-order derivative variation for the state.
 
void evalDpDm_diff (int n, Vec xbar, double Jbar)
 Computes derivative of second-order derivative variation term.
 
double evalEnergy (double time)
 Evaluates energy integral term.
 
void evalEnergy_diff (double time, double Jbar, Vec redgrad)
 Computes derivative of energy integral.
 

Protected Attributes

Vec stage
 
Vec stage_adj
 Intermediate stage vectors for forward and adjoint.
 
Vec rhs
 
Vec rhs_adj
 Right-hand side vectors for forward and adjoint.
 
KSP ksp
 PETSc's linear solver context for GMRES.
 
PC preconditioner
 Preconditioner for linear solver.
 
LinearSolverType linsolve_type
 Linear solver type (GMRES or NEUMANN)
 
int linsolve_maxiter
 Maximum number of linear solver iterations.
 
double linsolve_abstol
 Absolute tolerance for linear solver.
 
double linsolve_reltol
 Relative tolerance for linear solver.
 
int linsolve_iterstaken_avg
 Average number of linear solver iterations.
 
double linsolve_error_avg
 Average error of linear solver.
 
int linsolve_counter
 Counter for linear solve calls.
 
Vec tmp
 
Vec err
 Auxiliary vectors for Neumann iterations.
 
- Protected Attributes inherited from TimeStepper
Vec x
 Auxiliary vector for forward time stepping.
 
Vec xadj
 Auxiliary vector needed for adjoint (backward) time stepping.
 
Vec xprimal
 Auxiliary vector for backward time stepping.
 
std::vector< Vec > store_states
 Storage for primal states during forward evolution.
 
std::vector< Vec > dpdm_states
 Storage for states needed for second-order derivative penalty.
 
int mpirank_world
 MPI rank in global communicator.
 
int mpisize_petsc
 MPI size in Petsc communicator.
 
int mpirank_petsc
 MPI rank in Petsc communicator.
 
PetscInt localsize_u
 Size of local sub vector u or v in state x=[u,v].
 
PetscInt ilow
 First index of the local sub vector u,v.
 
PetscInt iupp
 Last index (+1) of the local sub vector u,v.
 
bool eval_leakage
 Flag to compute leakage integral term.
 
bool eval_energy
 Flag to compute energy integral term.
 
bool eval_dpdm
 Flag to compute second-order derivative integral term.
 
bool eval_weightedcost
 Flag to compute weighted cost integral term.
 
double leakage_integral
 Sums the integral over leakage.
 
double energy_integral
 Sums the energy term.
 
double dpdm_integral
 Sums second-order derivative variation value.
 
double weightedcost_integral
 Sums the integral over weighted cost function.
 
double weightedcost_width
 Width parameter for weighted cost function.
 

Additional Inherited Members

- Public Attributes inherited from TimeStepper
MasterEqmastereq
 Pointer to master equation solver.
 
int ntime
 Number of time steps.
 
double total_time
 Final evolution time.
 
double dt
 Time step size.
 
bool writeTrajectoryDataFiles
 Flag to determine whether or not trajectory data will be written to files during forward simulation *‍/.
 
Vec redgrad
 Reduced gradient vector for optimization.
 
OptimTargetoptim_target
 Pointer to optimization target specification.
 
Outputoutput
 Pointer to output handler.
 
bool storeFWD
 Flag to store primal states during forward evaluation.
 

Detailed Description

Implicit midpoint rule time integration scheme.

Second-order implicit method with symplectic properties. This is the default and recommended time integration, due to its stability and energy conservation properties.

Runge-Kutta tableau:

1/2 | 1/2
----------
| 1

Constructor & Destructor Documentation

◆ ImplMidpoint()

ImplMidpoint::ImplMidpoint ( MasterEq mastereq_,
int  ntime_,
double  total_time_,
LinearSolverType  linsolve_type_,
int  linsolve_maxiter_,
Output output_,
bool  storeFWD_ 
)

Constructor for implicit midpoint scheme.

Parameters
mastereq_Pointer to master equation solver
ntime_Number of time steps
total_time_Final evolution time
linsolve_type_Linear solver type (GMRES or NEUMANN)
linsolve_maxiter_Maximum linear solver iterations
output_Pointer to output handler
storeFWD_Flag to store forward states

◆ ~ImplMidpoint()

ImplMidpoint::~ImplMidpoint ( )

Member Function Documentation

◆ evolveBWD()

void ImplMidpoint::evolveBWD ( const double  tstart,
const double  tstop,
const Vec  x_stop,
Vec  x_adj,
Vec  grad,
bool  compute_gradient 
)
virtual

Evolves adjoint backward using implicit midpoint rule and adds to reduced gradient.

Parameters
tstartStart time (backward evolution)
tstopStop time (backward evolution)
x_stopState at stop time
x_adjAdjoint state vector
gradGradient vector to update
compute_gradientFlag to compute gradient

Reimplemented from TimeStepper.

Reimplemented in CompositionalImplMidpoint.

◆ evolveFWD()

void ImplMidpoint::evolveFWD ( const double  tstart,
const double  tstop,
Vec  x 
)
virtual

Evolves state forward using implicit midpoint rule.

Parameters
tstartStart time
tstopStop time
xState vector to evolve

Implements TimeStepper.

Reimplemented in CompositionalImplMidpoint.

◆ NeumannSolve()

int ImplMidpoint::NeumannSolve ( Mat  A,
Vec  b,
Vec  x,
double  alpha,
bool  transpose 
)

Solves (I - alpha*A) * x = b using Neumann iterations.

Parameters
AMatrix A
bRight-hand side vector
xSolution vector
alphaScaling parameter
transposeFlag to solve transposed system (I - alpha*A^T)*x = b
Returns
int Number of iterations taken

Member Data Documentation

◆ err

Vec ImplMidpoint::err
protected

Auxiliary vectors for Neumann iterations.

◆ ksp

KSP ImplMidpoint::ksp
protected

PETSc's linear solver context for GMRES.

◆ linsolve_abstol

double ImplMidpoint::linsolve_abstol
protected

Absolute tolerance for linear solver.

◆ linsolve_counter

int ImplMidpoint::linsolve_counter
protected

Counter for linear solve calls.

◆ linsolve_error_avg

double ImplMidpoint::linsolve_error_avg
protected

Average error of linear solver.

◆ linsolve_iterstaken_avg

int ImplMidpoint::linsolve_iterstaken_avg
protected

Average number of linear solver iterations.

◆ linsolve_maxiter

int ImplMidpoint::linsolve_maxiter
protected

Maximum number of linear solver iterations.

◆ linsolve_reltol

double ImplMidpoint::linsolve_reltol
protected

Relative tolerance for linear solver.

◆ linsolve_type

LinearSolverType ImplMidpoint::linsolve_type
protected

Linear solver type (GMRES or NEUMANN)

◆ preconditioner

PC ImplMidpoint::preconditioner
protected

Preconditioner for linear solver.

◆ rhs

Vec ImplMidpoint::rhs
protected

◆ rhs_adj

Vec ImplMidpoint::rhs_adj
protected

Right-hand side vectors for forward and adjoint.

◆ stage

Vec ImplMidpoint::stage
protected

◆ stage_adj

Vec ImplMidpoint::stage_adj
protected

Intermediate stage vectors for forward and adjoint.

◆ tmp

Vec ImplMidpoint::tmp
protected

The documentation for this class was generated from the following files: