MUQ  0.4.3
muq::Modeling::LinearSDE Class Reference

Defines a linear time invariant stochastic differential equation with Gaussian process noise. More...

#include <LinearSDE.h>

Detailed Description

Defines a linear time invariant stochastic differential equation with Gaussian process noise.

This class defines a LTI SDE of the form

\[ \frac{\partial f(t)}{\partial t} = F f(t) + L w(t), \]

where \(f(t)\) is the solution in \(\mathbb{R}^M\), \(F\) is an \(M\timesM\) matrix, \(L\) is an \(M\times N\) matrix, and \(w(t)\) is a white noise process with an \(N\timesN\) covariance matrix \(Q\).

It is possible to ignore the process noise by setting the matrix \(L\) to a nullptr in the constructor. This is used in the statespace representation of the constant Gaussian process covariance kernel.

Definition at line 27 of file LinearSDE.h.

Public Member Functions

template<typename EigenType1 , typename EigenType2 >
 LinearSDE (EigenType1 const &Fin, EigenType2 const &Lin, Eigen::MatrixXd const &Qin, boost::property_tree::ptree options)
 
template<typename Type1 , typename Type2 >
 LinearSDE (std::shared_ptr< Type1 > Fin, std::shared_ptr< Type2 > Lin, Eigen::MatrixXd const &Qin, boost::property_tree::ptree options)
 
 LinearSDE (std::shared_ptr< muq::Modeling::LinearOperator > Fin, std::shared_ptr< muq::Modeling::LinearOperator > Lin, Eigen::MatrixXd const &Qin, boost::property_tree::ptree options)
 
 LinearSDE (std::shared_ptr< muq::Modeling::LinearOperator > Fin, boost::property_tree::ptree options)
 
template<typename EigenRefVector >
Eigen::VectorXd EvolveState (EigenRefVector const &f0, double T) const
 
template<typename EigenRefVector1 , typename EigenRefVector2 >
void EvolveState (EigenRefVector1 const &f0, double T, EigenRefVector2 f) const
 
template<typename EigenRefVector , typename EigenRefMatrix >
std::pair< Eigen::VectorXd, Eigen::MatrixXd > EvolveDistribution (EigenRefVector const &mu0, EigenRefMatrix const &gamma0, double T) const
 
template<typename EigenRefVector1 , typename EigenRefMatrix1 , typename EigenRefVector2 , typename EigenRefMatrix2 >
void EvolveDistribution (EigenRefVector1 const &mu0, EigenRefMatrix1 const &gamma0, double T, EigenRefVector2 mu, EigenRefMatrix2 gamma) const
 
template<typename EigenRefVector , typename EigenRefMatrix >
std::pair< Eigen::VectorXd, Eigen::MatrixXd > EvolveDistribution (std::pair< EigenRefVector, EigenRefMatrix > const &muCov, double T) const
 
std::pair< Eigen::MatrixXd, Eigen::MatrixXd > Discretize (double deltaT)
 
std::shared_ptr< muq::Modeling::LinearOperatorGetF () const
 
std::shared_ptr< muq::Modeling::LinearOperatorGetL () const
 
Eigen::MatrixXd const & GetQ () const
 

Static Public Member Functions

static std::shared_ptr< LinearSDEConcatenate (std::vector< std::shared_ptr< LinearSDE >> const &sdes, boost::property_tree::ptree options=boost::property_tree::ptree())
 Combines the states of multiple SDEs into a single monolitch SDE. More...
 

Public Attributes

const int stateDim
 The dimension of the state variable \(f(t)\). More...
 

Constructor & Destructor Documentation

◆ LinearSDE() [1/4]

template<typename EigenType1 , typename EigenType2 >
muq::Modeling::LinearSDE::LinearSDE ( EigenType1 const &  Fin,
EigenType2 const &  Lin,
Eigen::MatrixXd const &  Qin,
boost::property_tree::ptree  options 
)
inline

Definition at line 33 of file LinearSDE.h.

◆ LinearSDE() [2/4]

template<typename Type1 , typename Type2 >
muq::Modeling::LinearSDE::LinearSDE ( std::shared_ptr< Type1 >  Fin,
std::shared_ptr< Type2 >  Lin,
Eigen::MatrixXd const &  Qin,
boost::property_tree::ptree  options 
)
inline

Definition at line 43 of file LinearSDE.h.

◆ LinearSDE() [3/4]

LinearSDE::LinearSDE ( std::shared_ptr< muq::Modeling::LinearOperator Fin,
std::shared_ptr< muq::Modeling::LinearOperator Lin,
Eigen::MatrixXd const &  Qin,
boost::property_tree::ptree  options 
)

Definition at line 15 of file LinearSDE.cpp.

References ExtractOptions(), F, L, LQLT, Q, sqrtQ, and nlohmann::to_string().

◆ LinearSDE() [4/4]

LinearSDE::LinearSDE ( std::shared_ptr< muq::Modeling::LinearOperator Fin,
boost::property_tree::ptree  options 
)

Definition at line 45 of file LinearSDE.cpp.

Member Function Documentation

◆ Concatenate()

std::shared_ptr< LinearSDE > LinearSDE::Concatenate ( std::vector< std::shared_ptr< LinearSDE >> const &  sdes,
boost::property_tree::ptree  options = boost::property_tree::ptree() 
)
static

Combines the states of multiple SDEs into a single monolitch SDE.

Consider \(N\) different stochastic differential equations defined by matrices \(F_i\), \(L_i\), and process covariances \(Q_i\). This function creates a new SDE defined by block diagonal matrices \(F\), \(L\), and \(Q\):

\[ F = \left[\begin{array}{cccc} F_1 & 0 & \cdots & \\ 0 & F_2 & 0 \\ \vdots & & \ddots & \\ 0 & \cdots & & F_N \end{array}\right] \]

\[ L = \left[\begin{array}{cccc} L_1 & 0 & \cdots & \\ 0 & L_2 & 0 \\ \vdots & & \ddots & \\ 0 & \cdots & & L_N \end{array}\right] \]

\[ Q = \left[\begin{array}{cccc} Q_1 & 0 & \cdots & \\ 0 & Q_2 & 0 \\ \vdots & & \ddots & \\ 0 & \cdots & & Q_N \end{array}\right] \]

Definition at line 180 of file LinearSDE.cpp.

References F, L, Q, and stateDim.

Referenced by muq::Approximation::SumKernel::GetStateSpace().

◆ Discretize()

std::pair< Eigen::MatrixXd, Eigen::MatrixXd > LinearSDE::Discretize ( double  deltaT)

Compute a matrix A and covariance Q such that \(x(t+\delta t) = A x(t) + q\) where \(q\) is a normal random variable with covariance \(Q\).

Definition at line 213 of file LinearSDE.cpp.

References dt, F, integratorType, LQLT, and stateDim.

◆ EvolveDistribution() [1/3]

template<typename EigenRefVector , typename EigenRefMatrix >
std::pair<Eigen::VectorXd, Eigen::MatrixXd> muq::Modeling::LinearSDE::EvolveDistribution ( EigenRefVector const &  mu0,
EigenRefMatrix const &  gamma0,
double  T 
) const
inline

Given the mean and covariance of the solution at time \(t\), compute the mean and covariance of the solution at time \(t+T\).

Definition at line 110 of file LinearSDE.h.

Referenced by EvolveDistribution().

◆ EvolveDistribution() [2/3]

template<typename EigenRefVector1 , typename EigenRefMatrix1 , typename EigenRefVector2 , typename EigenRefMatrix2 >
void muq::Modeling::LinearSDE::EvolveDistribution ( EigenRefVector1 const &  mu0,
EigenRefMatrix1 const &  gamma0,
double  T,
EigenRefVector2  mu,
EigenRefMatrix2  gamma 
) const
inline

Definition at line 122 of file LinearSDE.h.

References dt, F, integratorType, and LQLT.

◆ EvolveDistribution() [3/3]

template<typename EigenRefVector , typename EigenRefMatrix >
std::pair<Eigen::VectorXd, Eigen::MatrixXd> muq::Modeling::LinearSDE::EvolveDistribution ( std::pair< EigenRefVector, EigenRefMatrix > const &  muCov,
double  T 
) const
inline

Evolve the mean and covariance of the system using a std::pair to hold the distribution.

Definition at line 212 of file LinearSDE.h.

References EvolveDistribution().

◆ EvolveState() [1/2]

template<typename EigenRefVector >
Eigen::VectorXd muq::Modeling::LinearSDE::EvolveState ( EigenRefVector const &  f0,
double  T 
) const
inline

Given \(f(t)\), the state of the system at time \(t\), return a random realization of the state at time \(t+\delta t\).

Definition at line 62 of file LinearSDE.h.

◆ EvolveState() [2/2]

template<typename EigenRefVector1 , typename EigenRefVector2 >
void muq::Modeling::LinearSDE::EvolveState ( EigenRefVector1 const &  f0,
double  T,
EigenRefVector2  f 
) const
inline

Definition at line 71 of file LinearSDE.h.

References dt, F, muq::Utilities::RandomGenerator::GetNormal(), L, and sqrtQ.

◆ ExtractOptions()

void LinearSDE::ExtractOptions ( boost::property_tree::ptree  options)
protected

Definition at line 49 of file LinearSDE.cpp.

References dt, nlohmann::detail::dtoa_impl::e, and integratorType.

Referenced by LinearSDE().

◆ GetF()

std::shared_ptr<muq::Modeling::LinearOperator> muq::Modeling::LinearSDE::GetF ( ) const
inline

Definition at line 243 of file LinearSDE.h.

References F.

◆ GetL()

std::shared_ptr<muq::Modeling::LinearOperator> muq::Modeling::LinearSDE::GetL ( ) const
inline

Definition at line 244 of file LinearSDE.h.

References L.

◆ GetQ()

Eigen::MatrixXd const& muq::Modeling::LinearSDE::GetQ ( ) const
inline

Definition at line 245 of file LinearSDE.h.

References Q.

Member Data Documentation

◆ dt

double muq::Modeling::LinearSDE::dt
protected

Definition at line 258 of file LinearSDE.h.

Referenced by Discretize(), EvolveDistribution(), EvolveState(), and ExtractOptions().

◆ F

std::shared_ptr<muq::Modeling::LinearOperator> muq::Modeling::LinearSDE::F
protected

Definition at line 252 of file LinearSDE.h.

Referenced by Concatenate(), Discretize(), EvolveDistribution(), EvolveState(), GetF(), and LinearSDE().

◆ integratorType

std::string muq::Modeling::LinearSDE::integratorType
protected

Definition at line 259 of file LinearSDE.h.

Referenced by Discretize(), EvolveDistribution(), and ExtractOptions().

◆ L

std::shared_ptr<muq::Modeling::LinearOperator> muq::Modeling::LinearSDE::L
protected

Definition at line 253 of file LinearSDE.h.

Referenced by Concatenate(), EvolveState(), GetL(), and LinearSDE().

◆ LQLT

Eigen::MatrixXd muq::Modeling::LinearSDE::LQLT
protected

Definition at line 261 of file LinearSDE.h.

Referenced by Discretize(), EvolveDistribution(), and LinearSDE().

◆ Q

Eigen::MatrixXd muq::Modeling::LinearSDE::Q
protected

Definition at line 255 of file LinearSDE.h.

Referenced by Concatenate(), GetQ(), and LinearSDE().

◆ sqrtQ

Eigen::MatrixXd muq::Modeling::LinearSDE::sqrtQ
protected

Definition at line 256 of file LinearSDE.h.

Referenced by EvolveState(), and LinearSDE().

◆ stateDim

const int muq::Modeling::LinearSDE::stateDim

The dimension of the state variable \(f(t)\).

Definition at line 240 of file LinearSDE.h.

Referenced by Concatenate(), and Discretize().


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