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

Base class for time integration schemes to evolve the quantum dynamics. More...

#include <timestepper.hpp>

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

Public Member Functions

 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.
 
virtual void evolveFWD (const double tstart, const double tstop, Vec x)=0
 Evolves state forward by one time-step from tstart to tstop.
 
virtual void evolveBWD (const double tstart, const double tstop, const Vec x_stop, Vec x_adj, Vec grad, bool compute_gradient)
 Evolves adjoint state backward by one time-step and updates reduced gradient.
 

Public Attributes

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.
 

Protected Attributes

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.
 

Detailed Description

Base class for time integration schemes to evolve the quantum dynamics.

This abstract class provides the interface for time-stepping methods used to integrate the quantum evolution equations (Lindblad master equation or Schroedinger equation). It supports both forward and adjoint time integration, handles output of evolution data, and evaluates penalty integral terms.

Main functionality:

This class contains references to:

Constructor & Destructor Documentation

◆ TimeStepper() [1/2]

TimeStepper::TimeStepper ( )

◆ TimeStepper() [2/2]

TimeStepper::TimeStepper ( MasterEq mastereq_,
int  ntime_,
double  total_time_,
Output output_,
bool  storeFWD_ 
)

Constructor for time stepper.

Parameters
mastereq_Pointer to master equation solver
ntime_Number of time steps
total_time_Final evolution time
output_Pointer to output handler
storeFWD_Flag to store forward states

◆ ~TimeStepper()

TimeStepper::~TimeStepper ( )
virtual

Member Function Documentation

◆ evalDpDm()

double TimeStepper::evalDpDm ( Vec  x,
Vec  xm1,
Vec  xm2 
)

Evaluates second-order derivative variation for the state.

Parameters
xCurrent state vector
xm1State vector at previous time step
xm2State vector at two time steps ago
Returns
double Second-order variation value

◆ evalDpDm_diff()

void TimeStepper::evalDpDm_diff ( int  n,
Vec  xbar,
double  Jbar 
)

Computes derivative of second-order derivative variation term.

Parameters
nTime step index
xbarAdjoint state vector to update
JbarAdjoint of dpdm term

◆ evalEnergy()

double TimeStepper::evalEnergy ( double  time)

Evaluates energy integral term.

Parameters
timeCurrent time
Returns
double Energy value

◆ evalEnergy_diff()

void TimeStepper::evalEnergy_diff ( double  time,
double  Jbar,
Vec  redgrad 
)

Computes derivative of energy integral.

Parameters
timeCurrent time
JbarAdjoint of energy
redgradReduced gradient vector to update

◆ evalLeakage()

double TimeStepper::evalLeakage ( const Vec  x)

Evaluates leakage into guard levels.

Parameters
xCurrent state vector
Returns
double Leakage term value / ntime

◆ evalLeakage_diff()

void TimeStepper::evalLeakage_diff ( const Vec  x,
Vec  xbar,
double  Jbar 
)

Computes derivative of leakage term.

Parameters
xCurrent state vector
xbarAdjoint state vector to update
JbarAdjoint of leakage integral term

◆ evalWeightedCost()

double TimeStepper::evalWeightedCost ( double  time,
const Vec  x 
)

Evaluates the weighted cost function term.

Parameters
timeCurrent time
xCurrent state vector
Returns
double Weighted cost integral term value

◆ evalWeightedCost_diff()

void TimeStepper::evalWeightedCost_diff ( double  time,
const Vec  x,
Vec  xbar,
double  Jbar 
)

Computes derivative of weighted cost function term.

Parameters
timeCurrent time
xCurrent state vector
xbarAdjoint state vector to update
JbarAdjoint of weighted cost term

◆ evolveBWD()

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

Evolves adjoint state backward by one time-step and updates reduced gradient.

Abstract base-class implementation is empty. Derived classes that need backward time-stepping should implement this function.

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 in ExplEuler, ImplMidpoint, and CompositionalImplMidpoint.

◆ evolveFWD()

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

Evolves state forward by one time-step from tstart to tstop.

Pure virtual function to be implemented by the derived time-stepping classes.

Parameters
tstartStart time
tstopStop time
xState vector to evolve

Implemented in ExplEuler, ImplMidpoint, and CompositionalImplMidpoint.

◆ getDPDMIntegral()

double TimeStepper::getDPDMIntegral ( )
inline

◆ getEnergyIntegral()

double TimeStepper::getEnergyIntegral ( )
inline

◆ getLeakageIntegral()

double TimeStepper::getLeakageIntegral ( )
inline

◆ getState()

Vec TimeStepper::getState ( size_t  tindex)

Retrieves stored state at a specific time index.

Parameters
tindexTime step index
Returns
Vec State vector at the specified time

◆ getWeightedCostIntegral()

double TimeStepper::getWeightedCostIntegral ( )
inline

◆ setEvalDPDM()

void TimeStepper::setEvalDPDM ( bool  flag)
inline

◆ setEvalEnergy()

void TimeStepper::setEvalEnergy ( bool  flag)
inline

◆ setEvalLeakage()

void TimeStepper::setEvalLeakage ( bool  flag)
inline

◆ setEvalWeightedCost()

void TimeStepper::setEvalWeightedCost ( bool  flag,
double  width 
)
inline

◆ solveAdjointODE()

void TimeStepper::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.

This performs backward time-stepping to backpropagate an adjoint initial condition at final time (aka a terminal condtion) to time t=0, while accumulating the reduced gradient.

Parameters
rho_t0_barTerminal condition for adjoint state
finalstateFinal state from forward evolution
Jbar_leakageAdjoint of leakage integral term
Jbar_weightedcostAdjoint of weighted cost integral term
Jbar_dpdmAdjoint of second-order derivative variation
Jbar_energyAdjoint of energy integral term

◆ solveODE()

Vec TimeStepper::solveODE ( int  initid,
Vec  rho_t0 
)

Solves the ODE forward in time.

This performs the time-stepping to propagate an initial condition to the final time.

Parameters
initidInitial condition identifier
rho_t0Initial state vector
Returns
Vec Final state vector at time T

Member Data Documentation

◆ dpdm_integral

double TimeStepper::dpdm_integral
protected

Sums second-order derivative variation value.

◆ dpdm_states

std::vector<Vec> TimeStepper::dpdm_states
protected

Storage for states needed for second-order derivative penalty.

◆ dt

double TimeStepper::dt

Time step size.

◆ energy_integral

double TimeStepper::energy_integral
protected

Sums the energy term.

◆ eval_dpdm

bool TimeStepper::eval_dpdm
protected

Flag to compute second-order derivative integral term.

◆ eval_energy

bool TimeStepper::eval_energy
protected

Flag to compute energy integral term.

◆ eval_leakage

bool TimeStepper::eval_leakage
protected

Flag to compute leakage integral term.

◆ eval_weightedcost

bool TimeStepper::eval_weightedcost
protected

Flag to compute weighted cost integral term.

◆ ilow

PetscInt TimeStepper::ilow
protected

First index of the local sub vector u,v.

◆ iupp

PetscInt TimeStepper::iupp
protected

Last index (+1) of the local sub vector u,v.

◆ leakage_integral

double TimeStepper::leakage_integral
protected

Sums the integral over leakage.

◆ localsize_u

PetscInt TimeStepper::localsize_u
protected

Size of local sub vector u or v in state x=[u,v].

◆ mastereq

MasterEq* TimeStepper::mastereq

Pointer to master equation solver.

◆ mpirank_petsc

int TimeStepper::mpirank_petsc
protected

MPI rank in Petsc communicator.

◆ mpirank_world

int TimeStepper::mpirank_world
protected

MPI rank in global communicator.

◆ mpisize_petsc

int TimeStepper::mpisize_petsc
protected

MPI size in Petsc communicator.

◆ ntime

int TimeStepper::ntime

Number of time steps.

◆ optim_target

OptimTarget* TimeStepper::optim_target

Pointer to optimization target specification.

◆ output

Output* TimeStepper::output

Pointer to output handler.

◆ redgrad

Vec TimeStepper::redgrad

Reduced gradient vector for optimization.

◆ store_states

std::vector<Vec> TimeStepper::store_states
protected

Storage for primal states during forward evolution.

◆ storeFWD

bool TimeStepper::storeFWD

Flag to store primal states during forward evaluation.

◆ total_time

double TimeStepper::total_time

Final evolution time.

◆ weightedcost_integral

double TimeStepper::weightedcost_integral
protected

Sums the integral over weighted cost function.

◆ weightedcost_width

double TimeStepper::weightedcost_width
protected

Width parameter for weighted cost function.

◆ writeTrajectoryDataFiles

bool TimeStepper::writeTrajectoryDataFiles

Flag to determine whether or not trajectory data will be written to files during forward simulation *‍/.

◆ x

Vec TimeStepper::x
protected

Auxiliary vector for forward time stepping.

◆ xadj

Vec TimeStepper::xadj
protected

Auxiliary vector needed for adjoint (backward) time stepping.

◆ xprimal

Vec TimeStepper::xprimal
protected

Auxiliary vector for backward time stepping.


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