riccaticpp
Loading...
Searching...
No Matches
riccati Namespace Reference

Namespaces

namespace  internal
 

Classes

class  arena_alloc
 
struct  arena_allocator
 
class  arena_matrix
 
struct  dummy_allocator
 
struct  is_complex
 
class  SolverInfo
 

Typedefs

template<typename Scalar >
using matrix_t = Eigen::Matrix<Scalar, -1, -1>
 
template<typename Scalar >
using vector_t = Eigen::Matrix<Scalar, -1, 1>
 
template<typename Scalar >
using row_vector_t = Eigen::Matrix<Scalar, 1, -1>
 
template<typename Scalar >
using array2d_t = Eigen::Matrix<Scalar, -1, -1>
 
template<typename Scalar >
using array1d_t = Eigen::Matrix<Scalar, -1, 1>
 
template<typename Scalar >
using row_array1d_t = Eigen::Matrix<Scalar, 1, -1>
 
template<typename T >
using require_not_floating_point = std::enable_if_t<!std::is_floating_point<std::decay_t<T>>::value>
 
template<typename T >
using require_floating_point = std::enable_if_t<std::is_floating_point<std::decay_t<T>>::value>
 
template<typename T >
using require_floating_point_or_complex
 
template<typename T >
using require_not_floating_point_or_complex
 

Functions

template<typename T , typename Expr >
auto to_arena (arena_allocator< T, arena_alloc > &arena, const Expr &expr) noexcept
 
template<typename T , typename Expr >
auto eval (arena_allocator< T, arena_alloc > &arena, const Expr &expr) noexcept
 
template<typename Expr >
auto to_arena (dummy_allocator &arena, const Expr &expr) noexcept
 
template<typename T >
void print (const char *name, const arena_matrix< T > &x)
 
template<typename Mat >
RICCATI_ALWAYS_INLINE auto coeffs_to_cheby_nodes (Mat &&coeffs)
 Convert the Chebyshev coefficient representation of a set of polynomials P_j to their values at Chebyshev nodes of the second kind.
 
template<typename Mat >
RICCATI_ALWAYS_INLINE auto cheby_nodes_to_coeffs (Mat &&values)
 Convert a matrix of values of m polynomials evaluated at n+1 Chebyshev nodes of the second kind to their interpolating Chebyshev coefficients.
 
template<typename Mat >
RICCATI_ALWAYS_INLINE auto coeffs_and_cheby_nodes (Mat &&values)
 
template<typename Scalar , typename Integral >
RICCATI_ALWAYS_INLINE auto integration_matrix (Integral n)
 Constructs a Chebyshev integration matrix.
 
template<typename Scalar , typename Integral >
RICCATI_ALWAYS_INLINE auto quad_weights (Integral n)
 Calculates Clenshaw-Curtis quadrature weights.
 
template<typename Scalar , typename Integral >
RICCATI_ALWAYS_INLINE auto chebyshev (Integral n)
 Computes the Chebyshev differentiation matrix and Chebyshev nodes.
 
template<typename Vec1 , typename Vec2 , typename Allocator >
RICCATI_ALWAYS_INLINE auto interpolate (Vec1 &&s, Vec2 &&t, Allocator &&alloc)
 Creates an interpolation matrix from an array of source nodes to target nodes.
 
template<typename SolverInfo , typename Scalar , typename YScalar , typename Integral >
RICCATI_ALWAYS_INLINE auto spectral_chebyshev (SolverInfo &&info, Scalar x0, Scalar h, YScalar y0, YScalar dy0, Integral niter)
 Applies a spectral collocation method based on Chebyshev nodes for solving differential equations.
 
template<typename SolverInfo , typename Scalar , typename Vec >
auto osc_evolve (SolverInfo &&info, Scalar xi, Scalar xf, std::complex< Scalar > yi, std::complex< Scalar > dyi, Scalar eps, Scalar epsilon_h, Scalar init_stepsize, Vec &&x_eval, bool hard_stop=false)
 Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.
 
template<typename SolverInfo , typename Scalar , typename Vec >
auto nonosc_evolve (SolverInfo &&info, Scalar xi, Scalar xf, std::complex< Scalar > yi, std::complex< Scalar > dyi, Scalar eps, Scalar epsilon_h, Scalar init_stepsize, Vec &&x_eval, bool hard_stop=false)
 Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.
 
template<typename SolverInfo , typename Scalar , typename Vec >
auto evolve (SolverInfo &info, Scalar xi, Scalar xf, std::complex< Scalar > yi, std::complex< Scalar > dyi, Scalar eps, Scalar epsilon_h, Scalar init_stepsize, Vec &&x_eval, bool hard_stop=false)
 Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.
 
template<typename Expr >
RICCATI_ALWAYS_INLINE auto eval (dummy_allocator &alloc, Expr &&expr)
 
template<typename SolverInfo , typename Scalar , std::enable_if_t<!std::is_same< typename std::decay_t< SolverInfo >::funtype, pybind11::object >::value > * = nullptr>
auto gamma (SolverInfo &&info, const Scalar &x)
 
template<typename SolverInfo , typename Scalar , std::enable_if_t<!std::is_same< typename std::decay_t< SolverInfo >::funtype, pybind11::object >::value > * = nullptr>
auto omega (SolverInfo &&info, const Scalar &x)
 
template<typename Scalar , typename OmegaFun , typename Allocator , typename GammaFun , typename Integral >
auto make_solver (OmegaFun &&omega_fun, GammaFun &&gamma_fun, Allocator &&alloc, Integral nini, Integral nmax, Integral n, Integral p)
 Construct a new Solver Info object.
 
template<typename SolverInfo , typename Scalar , typename YScalar >
auto nonosc_step (SolverInfo &&info, Scalar x0, Scalar h, YScalar y0, YScalar dy0, Scalar epsres)
 Performs a single Chebyshev step with adaptive node count for solving differential equations.
 
template<bool DenseOut, typename SolverInfo , typename OmegaVec , typename GammaVec , typename Scalar , typename YScalar >
auto osc_step (SolverInfo &&info, OmegaVec &&omega_s, GammaVec &&gamma_s, Scalar x0, Scalar h, YScalar y0, YScalar dy0, Scalar epsres)
 Performs a single Riccati step for solving differential equations with oscillatory behavior.
 
template<typename SolverInfo , typename Scalar , typename Vec >
auto step (SolverInfo &info, Scalar xi, Scalar xf, std::complex< Scalar > yi, std::complex< Scalar > dyi, Scalar eps, Scalar epsilon_h, Scalar init_stepsize, Vec &&x_eval, bool hard_stop=false)
 Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.
 
template<typename SolverInfo , typename FloatingPoint >
FloatingPoint choose_nonosc_stepsize (SolverInfo &&info, FloatingPoint x0, FloatingPoint h, FloatingPoint epsilon_h)
 
template<typename SolverInfo , typename FloatingPoint >
auto choose_osc_stepsize (SolverInfo &&info, FloatingPoint x0, FloatingPoint h, FloatingPoint epsilon_h)
 Chooses an appropriate step size for the Riccati step based on the accuracy of Chebyshev interpolation of w(x) and g(x).
 
template<typename Scalar , typename Vector >
auto scale (Vector &&x, Scalar x0, Scalar h)
 Scales and shifts a vector of Chebyshev nodes.
 
template<typename T >
constexpr T pi ()
 
double eval (double x)
 
template<typename T >
std::complex< T > & eval (std::complex< T > &x)
 
template<typename T >
std::complex< T > eval (std::complex< T > &&x)
 
template<typename T >
auto eval (T &&x)
 
template<typename T , Eigen::Index R, Eigen::Index C>
Eigen::Matrix< T, R, C > eval (Eigen::Matrix< T, R, C > &&x)
 
template<typename T , Eigen::Index R, Eigen::Index C>
auto & eval (Eigen::Matrix< T, R, C > &x)
 
template<typename T , Eigen::Index R, Eigen::Index C>
const auto & eval (const Eigen::Matrix< T, R, C > &x)
 
template<typename T , Eigen::Index R, Eigen::Index C>
Eigen::Array< T, R, C > eval (Eigen::Array< T, R, C > &&x)
 
template<typename T , Eigen::Index R, Eigen::Index C>
auto & eval (Eigen::Array< T, R, C > &x)
 
template<typename T , Eigen::Index R, Eigen::Index C>
const auto & eval (const Eigen::Array< T, R, C > &x)
 
template<typename T , typename Scalar >
auto get_slice (T &&x_eval, Scalar start, Scalar end)
 
template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto sin (T x)
 
template<typename T , require_not_floating_point< T > * = nullptr>
auto sin (T &&x)
 
template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto cos (T x)
 
template<typename T , require_not_floating_point< T > * = nullptr>
auto cos (T &&x)
 
template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto sqrt (T x)
 
template<typename T , require_not_floating_point< T > * = nullptr>
auto sqrt (T &&x)
 
template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto square (T x)
 
template<typename T , require_not_floating_point< T > * = nullptr>
auto square (T &&x)
 
template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto array (T x)
 
template<typename T , require_not_floating_point< T > * = nullptr>
auto array (T &&x)
 
template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto matrix (T x)
 
template<typename T , require_not_floating_point< T > * = nullptr>
auto matrix (T &&x)
 
template<typename T , require_floating_point_or_complex< T > * = nullptr>
constexpr T zero_like (T x)
 
template<typename T , require_not_floating_point< T > * = nullptr>
auto zero_like (const T &x)
 
template<typename T1 , typename T2 , require_floating_point< T1 > * = nullptr>
auto pow (T1 x, T2 y)
 
template<typename T1 , typename T2 , require_not_floating_point< T1 > * = nullptr>
auto pow (T1 &&x, T2 y)
 
template<typename T , int R, int C>
void print (const char *name, const Eigen::Matrix< T, R, C > &x)
 
template<typename T , int R, int C>
void print (const char *name, const Eigen::Array< T, R, C > &x)
 
template<typename T , std::enable_if_t< std::is_floating_point< std::decay_t< T > >::value||std::is_integral< std::decay_t< T > >::value > * = nullptr>
void print (const char *name, T &&x)
 
template<typename T , std::enable_if_t< std::is_floating_point< std::decay_t< T > >::value||std::is_integral< std::decay_t< T > >::value > * = nullptr>
void print (const char *name, const std::complex< T > &x)
 

Variables

template<typename T >
constexpr Eigen::Index compile_size_v = std::decay_t<T>::RowsAtCompileTime * std::decay_t<T>::ColsAtCompileTime
 

Typedef Documentation

◆ array1d_t

template<typename Scalar >
using riccati::array1d_t = Eigen::Matrix<Scalar, -1, 1>

Definition at line 58 of file utils.hpp.

◆ array2d_t

template<typename Scalar >
using riccati::array2d_t = Eigen::Matrix<Scalar, -1, -1>

Definition at line 56 of file utils.hpp.

◆ matrix_t

template<typename Scalar >
using riccati::matrix_t = Eigen::Matrix<Scalar, -1, -1>

Definition at line 48 of file utils.hpp.

◆ require_floating_point

template<typename T >
using riccati::require_floating_point = std::enable_if_t<std::is_floating_point<std::decay_t<T>>::value>

Definition at line 140 of file utils.hpp.

◆ require_floating_point_or_complex

template<typename T >
using riccati::require_floating_point_or_complex
Initial value:
std::enable_if_t<std::is_floating_point<std::decay_t<T>>::value
|| is_complex<std::decay_t<T>>::value>

Definition at line 153 of file utils.hpp.

◆ require_not_floating_point

template<typename T >
using riccati::require_not_floating_point = std::enable_if_t<!std::is_floating_point<std::decay_t<T>>::value>

Definition at line 137 of file utils.hpp.

◆ require_not_floating_point_or_complex

Initial value:
std::enable_if_t<!std::is_floating_point<std::decay_t<T>>::value
&& !is_complex<std::decay_t<T>>::value>

Definition at line 157 of file utils.hpp.

◆ row_array1d_t

template<typename Scalar >
using riccati::row_array1d_t = Eigen::Matrix<Scalar, 1, -1>

Definition at line 60 of file utils.hpp.

◆ row_vector_t

template<typename Scalar >
using riccati::row_vector_t = Eigen::Matrix<Scalar, 1, -1>

Definition at line 53 of file utils.hpp.

◆ vector_t

template<typename Scalar >
using riccati::vector_t = Eigen::Matrix<Scalar, -1, 1>

Definition at line 50 of file utils.hpp.

Function Documentation

◆ array() [1/2]

template<typename T , require_not_floating_point< T > * = nullptr>
auto riccati::array ( T && x)
inline

Definition at line 206 of file utils.hpp.

◆ array() [2/2]

template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto riccati::array ( T x)
inline

Definition at line 201 of file utils.hpp.

◆ cheby_nodes_to_coeffs()

template<typename Mat >
RICCATI_ALWAYS_INLINE auto riccati::cheby_nodes_to_coeffs ( Mat && values)

Convert a matrix of values of m polynomials evaluated at n+1 Chebyshev nodes of the second kind to their interpolating Chebyshev coefficients.

This function computes the Chebyshev coefficients for a set of polynomials. The input is a matrix V, where each column contains the values of a polynomial at Chebyshev nodes. The output is a matrix C, where C(i, j) is the coefficient of the i-th Chebyshev polynomial for the j-th input polynomial. The relationship is given by:

\[ F_j(x) = \sum_{k=0}^{n} C_{kj}T_k(x) \]

which interpolates the values [V_{0j}, V_{1j}, ..., V_{nj}] for j = 0...(m-1).

This implementation is adapted from the vals2coeffs function in the chebfun package (https://github.com/chebfun/chebfun/blob/master/%40chebtech2/vals2coeffs.m).

Template Parameters
MatEigen matrix type
Parameters
valuesA matrix of size (n+1, m), where the (i, j)th element is the value of the j-th polynomial evaluated at the i-th Chebyshev node.
Returns
A matrix of size (n+1, m), where the (i, j)th element is the coefficient of the i-th Chebyshev polynomial for interpolating the j-th input polynomial.

Definition at line 105 of file chebyshev.hpp.

References riccati::internal::fft().

◆ chebyshev()

template<typename Scalar , typename Integral >
RICCATI_ALWAYS_INLINE auto riccati::chebyshev ( Integral n)

Computes the Chebyshev differentiation matrix and Chebyshev nodes.

This function calculates the Chebyshev differentiation matrix D of size (n+1, n+1) and n+1 Chebyshev nodes x for the standard 1D interval [-1, 1]. The differentiation matrix D can be used to approximate the derivative of a function sampled at the Chebyshev nodes. The nodes are computed according to the formula:

\[ x_p = \cos \left( \frac{\pi p}{n} \right), \quad p = 0, 1, \ldots, n. \]

Template Parameters
ScalarThe scalar type of the Chebyshev nodes and differentiation
IntegralThe integral type of the number of Chebyshev nodes
Parameters
nint - The number of Chebyshev nodes minus one.
Returns
std::pair<Eigen::MatrixXd, Eigen::VectorXd> - A pair consisting of:
  1. The differentiation matrix D of size (n+1, n+1).
  2. Eigen::VectorXd (real) - The vector of Chebyshev nodes x of size (n+1), ordered in descending order from 1 to -1.

Definition at line 278 of file chebyshev.hpp.

References array(), eval(), matrix(), and pi().

◆ choose_nonosc_stepsize()

template<typename SolverInfo , typename FloatingPoint >
FloatingPoint riccati::choose_nonosc_stepsize ( SolverInfo && info,
FloatingPoint x0,
FloatingPoint h,
FloatingPoint epsilon_h )
inline

Chooses the stepsize for spectral Chebyshev steps, based on the variation of 1/w, the approximate timescale over which the solution changes. If over the suggested interval h 1/w changes by a fraction of

\[\pm epsilon_h\]

or more, the interval is halved, otherwise it's accepted.

Template Parameters
SolverInfoA riccati solver like object
FloatingPointA floating point
Parameters
infoSolverinfo object which is used to retrieve Solverinfo.xp, the (p+1) Chebyshev nodes used for interpolation to determine the stepsize.
x0Current value of the independent variable.
hInitial estimate of the stepsize.
epsilon_hTolerance parameter defining how much 1/w(x) is allowed to change over the course of the step.
Returns
Refined stepsize over which 1/w(x) does not change by more than epsilon_h/w(x).

Definition at line 27 of file stepsize.hpp.

References choose_nonosc_stepsize(), omega(), and scale().

◆ choose_osc_stepsize()

template<typename SolverInfo , typename FloatingPoint >
auto riccati::choose_osc_stepsize ( SolverInfo && info,
FloatingPoint x0,
FloatingPoint h,
FloatingPoint epsilon_h )
inline

Chooses an appropriate step size for the Riccati step based on the accuracy of Chebyshev interpolation of w(x) and g(x).

This function determines an optimal step size h over which the functions w(x) and g(x) can be represented with sufficient accuracy by evaluating their values at p+1 Chebyshev nodes. It performs interpolation to p points halfway between these nodes and compares the interpolated values with the actual values of w(x) and g(x). If the largest relative error in w or g exceeds the tolerance epsh, the step size h is reduced. This process ensures that the Chebyshev interpolation of w(x) and g(x) over the step [x0, x0+h] has a relative error no larger than epsh.

Parameters
infoSolverInfo object - Object containing pre-computed information and methods for evaluating functions w(x) and g(x), as well as interpolation matrices and node positions.
x0float - The current value of the independent variable.
hfloat - The initial estimate of the step size.
epsilon_hfloat - Tolerance parameter defining the maximum relative error allowed in the Chebyshev interpolation of w(x) and g(x) over the proposed step.
Returns
float - The refined step size over which the Chebyshev interpolation of w(x) and g(x) satisfies the relative error tolerance epsh.

Definition at line 64 of file stepsize.hpp.

References choose_osc_stepsize(), eval(), gamma(), omega(), and scale().

◆ coeffs_and_cheby_nodes()

template<typename Mat >
RICCATI_ALWAYS_INLINE auto riccati::coeffs_and_cheby_nodes ( Mat && values)

Returns both the Chebyshev coefficients and the Chebyshev nodes of the second kind for a set of polynomials.

Template Parameters
MatEigen matrix type
Parameters
valuesA matrix of size (n+1, m), where the (i, j)th element is the value of the j-th polynomial evaluated at the i-th Chebyshev node.

Definition at line 132 of file chebyshev.hpp.

References eval(), and riccati::internal::fft().

◆ coeffs_to_cheby_nodes()

template<typename Mat >
RICCATI_ALWAYS_INLINE auto riccati::coeffs_to_cheby_nodes ( Mat && coeffs)

Convert the Chebyshev coefficient representation of a set of polynomials P_j to their values at Chebyshev nodes of the second kind.

This function computes the values of a set of polynomials at Chebyshev nodes of the second kind. The input is a matrix of coefficients C, where each column represents a polynomial. The output is a matrix V, where V(i,j) is the value of the j-th polynomial at the i-th Chebyshev node. The relationship is given by:

\[ V_{ij} = P_j(x_i) = \sum_{k=0}^{n} C_{kj}T_k(x_i) \]

This implementation is adapted from the coeff2vals function in the chebfun package (https://github.com/chebfun/chebfun/blob/master/%40chebtech2/coeffs2vals.m).

Template Parameters
MatEigen matrix type
Parameters
coeffsA matrix of size (n+1, m), where the (i, j)th element represents the projection of the j-th input polynomial onto the i-th Chebyshev polynomial.
Returns
A matrix of size (n+1, m), where the (i, j)th element represents the j-th input polynomial evaluated at the i-th Chebyshev node.

Definition at line 58 of file chebyshev.hpp.

References eval(), and riccati::internal::fft().

◆ cos() [1/2]

template<typename T , require_not_floating_point< T > * = nullptr>
auto riccati::cos ( T && x)
inline

Definition at line 176 of file utils.hpp.

◆ cos() [2/2]

template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto riccati::cos ( T x)
inline

Definition at line 171 of file utils.hpp.

◆ eval() [1/12]

template<typename T , typename Expr >
auto riccati::eval ( arena_allocator< T, arena_alloc > & arena,
const Expr & expr )
inlinenoexcept

Definition at line 161 of file arena_matrix.hpp.

◆ eval() [2/12]

template<typename T , Eigen::Index R, Eigen::Index C>
const auto & riccati::eval ( const Eigen::Array< T, R, C > & x)
inline

Definition at line 108 of file utils.hpp.

◆ eval() [3/12]

template<typename T , Eigen::Index R, Eigen::Index C>
const auto & riccati::eval ( const Eigen::Matrix< T, R, C > & x)
inline

Definition at line 93 of file utils.hpp.

◆ eval() [4/12]

double riccati::eval ( double x)
inline

Definition at line 67 of file utils.hpp.

◆ eval() [5/12]

template<typename Expr >
RICCATI_ALWAYS_INLINE auto riccati::eval ( dummy_allocator & alloc,
Expr && expr )

Definition at line 396 of file memory.hpp.

References eval().

◆ eval() [6/12]

template<typename T , Eigen::Index R, Eigen::Index C>
Eigen::Array< T, R, C > riccati::eval ( Eigen::Array< T, R, C > && x)
inline

Definition at line 98 of file utils.hpp.

◆ eval() [7/12]

template<typename T , Eigen::Index R, Eigen::Index C>
auto & riccati::eval ( Eigen::Array< T, R, C > & x)
inline

Definition at line 103 of file utils.hpp.

◆ eval() [8/12]

template<typename T , Eigen::Index R, Eigen::Index C>
Eigen::Matrix< T, R, C > riccati::eval ( Eigen::Matrix< T, R, C > && x)
inline

Definition at line 83 of file utils.hpp.

◆ eval() [9/12]

template<typename T , Eigen::Index R, Eigen::Index C>
auto & riccati::eval ( Eigen::Matrix< T, R, C > & x)
inline

Definition at line 88 of file utils.hpp.

◆ eval() [10/12]

template<typename T >
std::complex< T > riccati::eval ( std::complex< T > && x)
inline

Definition at line 73 of file utils.hpp.

◆ eval() [11/12]

template<typename T >
std::complex< T > & riccati::eval ( std::complex< T > & x)
inline

Definition at line 69 of file utils.hpp.

◆ eval() [12/12]

template<typename T >
auto riccati::eval ( T && x)
inline

Definition at line 78 of file utils.hpp.

◆ evolve()

template<typename SolverInfo , typename Scalar , typename Vec >
auto riccati::evolve ( SolverInfo & info,
Scalar xi,
Scalar xf,
std::complex< Scalar > yi,
std::complex< Scalar > dyi,
Scalar eps,
Scalar epsilon_h,
Scalar init_stepsize,
Vec && x_eval,
bool hard_stop = false )
inline

Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.

This function solves the differential equation on the interval (xi, xf), starting from the initial conditions y(xi) = yi and y'(xi) = dyi. It keeps the residual of the ODE below eps, and returns an interpolated solution (dense output) at the points specified in x_eval.

Template Parameters
SolverInfoType of the solver info object containing differentiation matrices, etc.
ScalarNumeric scalar type, typically float or double.
VecType of the vector for dense output values, should match Scalar type.
SolverInfoType of the solver info object containing differentiation matrices, etc.
ScalarNumeric scalar type, typically float or double.
VecType of the vector for dense output values.
Parameters
[in]infoSolverInfo object containing necessary information for the solver.
[in]xiStarting value of the independent variable.
[in]xfEnding value of the independent variable.
[in]yiInitial value of the dependent variable at xi.
[in]dyiInitial derivative of the dependent variable at xi.
[in]epsRelative tolerance for the local error of both Riccati and Chebyshev type steps.
[in]epsilon_hRelative tolerance for choosing the stepsize of Riccati steps.
[in]init_stepsizeinitial stepsize for the integration
[in]x_evalList of x-values where the solution is to be interpolated (dense output) and returned.
[in]hard_stopIf true, forces the solver to have a potentially smaller last stepsize to stop exactly at xf.
Returns
A tuple containing multiple elements representing the results of the ODE solving process: 0. std::vector<Scalar>: A vector containing the x-values at which the solution was evaluated or interpolated. These values correspond to the points in the interval [xi, xf] and include the points specified in x_eval if dense output was requested.
  1. std::vector<std::complex<Scalar>>: A vector of complex numbers representing the solution y(x) of the differential equation at each x-value from the corresponding vector of x-values.
  2. std::vector<std::complex<Scalar>>: A vector of complex numbers representing the derivative of the solution, y'(x), at each x-value from the corresponding vector of x-values.
  3. std::vector<int>: A vector indicating the success status of the solver at each step. Each element corresponds to a step in the solver process, where a value of 1 indicates success, and 0 indicates failure.
  4. std::vector<int>: A vector indicating the type of step taken at each point in the solution process. Each element corresponds to a step in the solver process, where a value of 1 indicates an oscillatory step and 0 indicates a non-oscillatory step.
  5. std::vector<Scalar>: A vector containing the phase angle at each step of the solution process if relevant. This is especially applicable for oscillatory solutions and may not be used for all types of differential equations.
  6. Eigen::Matrix<std::complex<Scalar>, -1, 1>: A vector containing the interpolated solution at the specified x_eval The function returns these vectors encapsulated in a standard tuple, providing comprehensive information about the solution process, including where the solution was evaluated, the values and derivatives of the solution at those points, success status of each step, type of each step, and phase angles where applicable.

Definition at line 319 of file evolve.hpp.

References riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::alloc_, riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::chebyshev_, choose_nonosc_stepsize(), choose_osc_stepsize(), compile_size_v, riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::Dn(), eval(), gamma(), get_slice(), interpolate(), matrix(), nonosc_step(), omega(), osc_step(), pi(), scale(), riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::xn(), and riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::xp().

◆ gamma()

template<typename SolverInfo , typename Scalar , std::enable_if_t<!std::is_same< typename std::decay_t< SolverInfo >::funtype, pybind11::object >::value > * = nullptr>
auto riccati::gamma ( SolverInfo && info,
const Scalar & x )
inline

Definition at line 25 of file solver.hpp.

◆ get_slice()

template<typename T , typename Scalar >
auto riccati::get_slice ( T && x_eval,
Scalar start,
Scalar end )

Definition at line 113 of file utils.hpp.

◆ integration_matrix()

template<typename Scalar , typename Integral >
RICCATI_ALWAYS_INLINE auto riccati::integration_matrix ( Integral n)

Constructs a Chebyshev integration matrix.

This function computes the Chebyshev integration matrix, which maps the values of a function at n Chebyshev nodes of the second kind (ordered from +1 to -1) to the values of the integral of the interpolating polynomial at those nodes. The integral is computed on the interval defined by the Chebyshev nodes, with the last value of the integral (at the start of the interval) being set to zero. This implementation is adapted from the cumsummat function in the chebfun package (https://github.com/chebfun/chebfun/blob/master/%40chebcolloc2/chebcolloc2.m).

Template Parameters
ScalarThe scalar type of the Chebyshev nodes and integration matrix
IntegralThe integral type of the number of Chebyshev nodes
Parameters
nNumber of Chebyshev nodes the integrand is evaluated at. The nodes are ordered from +1 to -1.
Returns
Integration matrix of size (n, n). This matrix maps the values of the integrand at the n Chebyshev nodes to the values of the definite integral on the interval, up to each of the Chebyshev nodes (the last value being zero by definition).

Definition at line 178 of file chebyshev.hpp.

References array(), coeffs_and_cheby_nodes(), and eval().

◆ interpolate()

template<typename Vec1 , typename Vec2 , typename Allocator >
RICCATI_ALWAYS_INLINE auto riccati::interpolate ( Vec1 && s,
Vec2 && t,
Allocator && alloc )

Creates an interpolation matrix from an array of source nodes to target nodes.

This function constructs an interpolation matrix that maps function values known at source nodes s to estimated values at target nodes t. The computation is based on the Vandermonde matrix approach and the resulting matrix L applies the interpolation. The method is adapted from the implementation provided here (https://github.com/ahbarnett/BIE3D/blob/master/utils/interpmat_1d.m).

Template Parameters
Vec1An Eigen vector type
Vec2An Eigen vector type
AllocatorAn allocator type
Parameters
sA vector specifying the source nodes, at which the function values are known.
tA vector specifying the target nodes, at which the function values are to be interpolated.
allocAn allocator for the Eigen objects.
Returns
The interpolation matrix L. If s has size p and t has size q, then L has size (q, p). L takes function values at source points s and yields the function evaluated at target points t.

Definition at line 332 of file chebyshev.hpp.

References eval().

◆ make_solver()

template<typename Scalar , typename OmegaFun , typename Allocator , typename GammaFun , typename Integral >
auto riccati::make_solver ( OmegaFun && omega_fun,
GammaFun && gamma_fun,
Allocator && alloc,
Integral nini,
Integral nmax,
Integral n,
Integral p )
inline

Construct a new Solver Info object.

Template Parameters
ScalarA scalar type used for calculations
OmegaFunType of the frequency function. Must be able to take in and return scalars and vectors.
GammaFunType of the friction function. Must be able to take in and return scalars and vectors.
IntegralAn integral type
Parameters
omega_funFrequency function
gamma_funFriction function
niniMinimum number of Chebyshev nodes to use inside Chebyshev collocation steps.
nmaxMaximum number of Chebyshev nodes to use inside Chebyshev collocation steps.
n(Number of Chebyshev nodes - 1) to use inside Chebyshev collocation steps.
p(Number of Chebyshev nodes - 1) to use for computing Riccati steps.

Definition at line 274 of file solver.hpp.

◆ matrix() [1/2]

template<typename T , require_not_floating_point< T > * = nullptr>
auto riccati::matrix ( T && x)
inline

Definition at line 216 of file utils.hpp.

◆ matrix() [2/2]

template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto riccati::matrix ( T x)
inline

Definition at line 211 of file utils.hpp.

◆ nonosc_evolve()

template<typename SolverInfo , typename Scalar , typename Vec >
auto riccati::nonosc_evolve ( SolverInfo && info,
Scalar xi,
Scalar xf,
std::complex< Scalar > yi,
std::complex< Scalar > dyi,
Scalar eps,
Scalar epsilon_h,
Scalar init_stepsize,
Vec && x_eval,
bool hard_stop = false )
inline

Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.

This function solves the differential equation on the interval (xi, xf), starting from the initial conditions y(xi) = yi and y'(xi) = dyi. It keeps the residual of the ODE below eps, and returns an interpolated solution (dense output) at the points specified in x_eval.

Template Parameters
SolverInfoType of the solver info object containing differentiation matrices, etc.
ScalarNumeric scalar type, typically float or double.
VecType of the vector for dense output values, should match Scalar type.
VecType of the vector for dense output values.
Parameters
[in]infoSolverInfo object containing necessary information for the solver.
[in]xiStarting value of the independent variable.
[in]xfEnding value of the independent variable.
[in]yiInitial value of the dependent variable at xi.
[in]dyiInitial derivative of the dependent variable at xi.
[in]epsRelative tolerance for the local error of both Riccati and Chebyshev type steps.
[in]epsilon_hRelative tolerance for choosing the stepsize of Riccati steps.
[in]init_stepsizeStepsize for the integration
[in]x_evalList of x-values where the solution is to be interpolated (dense output) and returned.
[in]hard_stopIf true, forces the solver to have a potentially smaller last stepsize to stop exactly at xf.
Returns
A tuple containing the following elements:
  1. bool: A boolean flag indicating whether the step was successfully taken.
  2. Scalar: The next x-value after the integration step, indicating the new position in the integration domain.
  3. Scalar: The suggested next stepsize for further integration steps based on the current step's data.
  4. std::tuple: The result from the riccati::nonosc_step function which includes various internal states and calculations specific to the current step.
  5. Eigen::Matrix<std::complex<Scalar>, -1, 1>: A vector containing the interpolated solution at the specified x_eval points. This represents the dense output for the current step.
  6. Eigen::Index: The starting index in the x_eval array corresponding to the first point in the current step's dense output.
  7. Eigen::Index: The number of points in the x_eval array that are covered by the current step's dense output.
Note
Tuple element 4 vector contains the interpolated values of the differential equation's solution at the points specified by the x_eval input parameter. These values are calculated using the dense output methodology and are meant for high-accuracy interpolation between the standard discrete steps of the solver.

Definition at line 199 of file evolve.hpp.

References choose_nonosc_stepsize(), compile_size_v, get_slice(), interpolate(), matrix(), nonosc_step(), and omega().

◆ nonosc_step()

template<typename SolverInfo , typename Scalar , typename YScalar >
auto riccati::nonosc_step ( SolverInfo && info,
Scalar x0,
Scalar h,
YScalar y0,
YScalar dy0,
Scalar epsres )
inline

Performs a single Chebyshev step with adaptive node count for solving differential equations.

This function advances the solution of a differential equation from x = x0 by a step of size h, starting from the initial conditions y(x0) = y0 and ‘y’(x0) = dy0. It employs a Chebyshev spectral method with an adaptive number of nodes, starting withinfo.nininodes and doubling the count in each iteration until the relative accuracyepsresis achieved or the number of nodes exceedsinfo.nmax. The relative error is assessed by comparing the predicted values of the dependent variable at the end of the step for the current and previous iterations. If the desired accuracy isn't met with the maximum number of nodes, step sizehmay need to be reduced,info.nmax` increased, or a different numerical method considered.

Parameters
infoSolverInfo object - Object containing pre-computed information for the solver, like differentiation matrices and Chebyshev nodes, and methods for evaluating functions w(x) and g(x) over the interval [x0, x0+h].
x0float - The starting value of the independent variable.
hfloat - Step size for the spectral method.
y0complex - Initial value of the dependent variable at x0.
dy0complex - Initial derivative of the dependent variable at x0.
epsresfloat - Tolerance for the relative accuracy of the Chebyshev step.
Returns
std::tuple<std::complex<double>, std::complex<double>, float, int> - A tuple containing:
  1. std::complex<double> - Value of the dependent variable at the end of the step, at x = x0 + h.
  2. std::complex<double> - Value of the derivative of the dependent variable at the end of the step, at x = x0 + h.
  3. float - The (absolute) value of the relative difference of the dependent variable at the end of the step as predicted in the last and the previous iteration.
  4. int - Flag indicating success (1) if the asymptotic series has reached the desired epsres residual, 0 otherwise.

Definition at line 50 of file step.hpp.

References spectral_chebyshev().

◆ omega()

template<typename SolverInfo , typename Scalar , std::enable_if_t<!std::is_same< typename std::decay_t< SolverInfo >::funtype, pybind11::object >::value > * = nullptr>
auto riccati::omega ( SolverInfo && info,
const Scalar & x )
inline

Definition at line 33 of file solver.hpp.

◆ osc_evolve()

template<typename SolverInfo , typename Scalar , typename Vec >
auto riccati::osc_evolve ( SolverInfo && info,
Scalar xi,
Scalar xf,
std::complex< Scalar > yi,
std::complex< Scalar > dyi,
Scalar eps,
Scalar epsilon_h,
Scalar init_stepsize,
Vec && x_eval,
bool hard_stop = false )
inline

Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.

This function solves the differential equation on the interval (xi, xf), starting from the initial conditions y(xi) = yi and y'(xi) = dyi. It keeps the residual of the ODE below eps, and returns an interpolated solution (dense output) at the points specified in x_eval.

Template Parameters
SolverInfoType of the solver info object containing differentiation matrices, etc.
ScalarNumeric scalar type, typically float or double.
VecType of the vector for dense output values, should match Scalar type.
VecType of the vector for dense output values.
AllocatorType of the allocator for the arena memory pool.
Parameters
[in]infoSolverInfo object containing necessary information for the solver.
[in]xiStarting value of the independent variable.
[in]xfEnding value of the independent variable.
[in]yiInitial value of the dependent variable at xi.
[in]dyiInitial derivative of the dependent variable at xi.
[in]epsRelative tolerance for the local error of both Riccati and Chebyshev type steps.
[in]epsilon_hRelative tolerance for choosing the stepsize of Riccati steps.
[in]init_stepsizeStepsize for the integration
[in]x_evalList of x-values where the solution is to be interpolated (dense output) and returned.
[in]hard_stopIf true, forces the solver to have a potentially smaller last stepsize to stop exactly at xf.
Returns
A tuple containing the following elements:
  1. bool: A boolean flag indicating whether the step was successfully taken.
  2. Scalar: The next x-value after the integration step, indicating the new position in the integration domain.
  3. Scalar: The suggested next stepsize for further integration steps based on the current step's data.
  4. std::tuple: The result from the riccati::osc_step function which includes various internal states and calculations specific to the current step.
  5. Eigen::Matrix<std::complex<Scalar>, -1, 1>: A vector containing the interpolated solution at the specified x_eval points. This represents the dense output for the current step.
  6. Eigen::Index: The starting index in the x_eval array corresponding to the first point in the current step's dense output.
  7. Eigen::Index: The number of points in the x_eval array that are covered by the current step's dense output.
Note
Tuple element 4 vector contains the interpolated values of the differential equation's solution at the points specified by the x_eval input parameter. These values are calculated using the dense output methodology and are meant for high-accuracy interpolation between the standard discrete steps of the solver.

Definition at line 69 of file evolve.hpp.

References choose_osc_stepsize(), compile_size_v, eval(), gamma(), get_slice(), interpolate(), omega(), osc_step(), and scale().

◆ osc_step()

template<bool DenseOut, typename SolverInfo , typename OmegaVec , typename GammaVec , typename Scalar , typename YScalar >
auto riccati::osc_step ( SolverInfo && info,
OmegaVec && omega_s,
GammaVec && gamma_s,
Scalar x0,
Scalar h,
YScalar y0,
YScalar dy0,
Scalar epsres )
inline

Performs a single Riccati step for solving differential equations with oscillatory behavior.

This function advances the solution from x0 by h, starting from the initial conditions y(x0) = y0 and ‘y’(x0) = dy0, using an asymptotic expansion approach tailored for Riccati-type differential equations. It iteratively increases the order of the asymptotic series used for the Riccati equation until a residual ofepsresis reached or the residual stops decreasing, indicating that the asymptotic series cannot approximate the solution with the required accuracy over the given interval. In such cases, it is recommended to reduce the step sizeh` or consider an alternative approximation method. The function also computes the total phase change of the dependent variable over the step.

Template Parameters
SolverInfoA riccati solver like object
OmegaVecAn Eigen vector
GammaVecAn Eigen vector
ScalarA scalar type for the x values
YScalarA scalar type for the y values
Parameters
infoSolverInfo object - Object containing pre-computed information for the solver, like differentiation matrices and methods for evaluating functions w(x) and g(x) over the interval [x0, x0+h].
omega_sVector of the frequency function w(x) evaluated at the Chebyshev nodes over the interval [x0, x0+h].
gamma_sVector of the friction function g(x) evaluated at the Chebyshev nodes over the interval [x0, x0+h].
x0float - The starting value of the independent variable.
hfloat - Step size for the Riccati step.
y0complex - Initial value of the dependent variable at x0.
dy0complex - Initial derivative of the dependent variable at x0.
epsresfloat - Tolerance for the relative accuracy of the Riccati step.
Returns
std::tuple<std::complex<double>, std::complex<double>, float, int, std::complex<double>> - A tuple containing:
  1. std::complex<double> - Value of the dependent variable at the end of the step, at x = x0 + h.
  2. std::complex<double> - Value of the derivative of the dependent variable at the end of the step, at x = x0 + h.
  3. float - Maximum value of the residual (after the final iteration of the asymptotic approximation) over the Chebyshev nodes across the interval.
  4. int - Flag indicating success (1) if the asymptotic series has reached the desired epsres residual, 0 otherwise.
  5. std::complex<double> - Total phase change (not mod 2π) of the dependent variable over the step.
Warning
This function relies on info.wn, info.gn being set correctly for a step of size h. If solve() is calling this function, that is taken care of automatically, but it needs to be done manually otherwise.

Definition at line 138 of file step.hpp.

References eval().

◆ pi()

template<typename T >
T riccati::pi ( )
inlineconstexpr

Definition at line 63 of file utils.hpp.

◆ pow() [1/2]

template<typename T1 , typename T2 , require_not_floating_point< T1 > * = nullptr>
auto riccati::pow ( T1 && x,
T2 y )
inline

Definition at line 236 of file utils.hpp.

◆ pow() [2/2]

template<typename T1 , typename T2 , require_floating_point< T1 > * = nullptr>
auto riccati::pow ( T1 x,
T2 y )
inline

Definition at line 231 of file utils.hpp.

◆ print() [1/5]

template<typename T >
void riccati::print ( const char * name,
const arena_matrix< T > & x )
inline

Definition at line 172 of file arena_matrix.hpp.

◆ print() [2/5]

template<typename T , int R, int C>
void riccati::print ( const char * name,
const Eigen::Array< T, R, C > & x )
inline

Definition at line 249 of file utils.hpp.

◆ print() [3/5]

template<typename T , int R, int C>
void riccati::print ( const char * name,
const Eigen::Matrix< T, R, C > & x )
inline

Definition at line 241 of file utils.hpp.

◆ print() [4/5]

template<typename T , std::enable_if_t< std::is_floating_point< std::decay_t< T > >::value||std::is_integral< std::decay_t< T > >::value > * = nullptr>
void riccati::print ( const char * name,
const std::complex< T > & x )
inline

Definition at line 270 of file utils.hpp.

◆ print() [5/5]

template<typename T , std::enable_if_t< std::is_floating_point< std::decay_t< T > >::value||std::is_integral< std::decay_t< T > >::value > * = nullptr>
void riccati::print ( const char * name,
T && x )
inline

Definition at line 260 of file utils.hpp.

◆ quad_weights()

template<typename Scalar , typename Integral >
RICCATI_ALWAYS_INLINE auto riccati::quad_weights ( Integral n)

Calculates Clenshaw-Curtis quadrature weights.

This function computes the Clenshaw-Curtis quadrature weights, which map function evaluations at n+1 Chebyshev nodes of the second kind (ordered from +1 to -1) to the value of the definite integral of the interpolating function on the same interval. The method is based on the Clenshaw-Curtis quadrature formula, as described in Trefethen's "Spectral Methods in MATLAB" (Chapter 12, clencurt.m).

Template Parameters
ScalarThe scalar type of the quadrature weights
IntegralThe integral type of the number of Chebyshev nodes
Parameters
nThe number of Chebyshev nodes minus one, for which the quadrature weights are to be computed.
Returns
A vector of size (n+1), containing the quadrature weights.
Note
See the below for more information Trefethen, Lloyd N. Spectral methods in MATLAB. Society for industrial and applied mathematics, 2000.

Definition at line 225 of file chebyshev.hpp.

References array(), and pi().

◆ scale()

template<typename Scalar , typename Vector >
auto riccati::scale ( Vector && x,
Scalar x0,
Scalar h )
inline

Scales and shifts a vector of Chebyshev nodes.

This function takes a vector of Chebyshev nodes and scales it to fit a specific interval. Each element x(i) of the input vector is transformed to fit into the new interval starting at x0 with a width of h, effectively mapping the standard Chebyshev interval [-1, 1] to [x0, x0 + h]. The transformation used is x(i) -> x0 + h/2 + h/2 * x(i). This is commonly used to adjust Chebyshev nodes or the results of spectral methods to the specific interval of interest in numerical computations.

Template Parameters
ScalarThe scalar type, typically a floating-point type like double or float.
VectorThe Eigen vector type, typically Eigen::VectorXd or similar.
Parameters
xVector (forwarded reference) - The input vector, typically containing Chebyshev nodes or similar quantities.
x0Scalar - The start of the interval to which the nodes should be scaled.
hScalar - The width of the interval to which the nodes should be scaled.
Returns
Returns a new Eigen vector of the same type as x, with each element scaled and shifted to the new interval [x0, x0 + h].

Definition at line 40 of file utils.hpp.

◆ sin() [1/2]

template<typename T , require_not_floating_point< T > * = nullptr>
auto riccati::sin ( T && x)
inline

Definition at line 166 of file utils.hpp.

◆ sin() [2/2]

template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto riccati::sin ( T x)
inline

Definition at line 161 of file utils.hpp.

◆ spectral_chebyshev()

template<typename SolverInfo , typename Scalar , typename YScalar , typename Integral >
RICCATI_ALWAYS_INLINE auto riccati::spectral_chebyshev ( SolverInfo && info,
Scalar x0,
Scalar h,
YScalar y0,
YScalar dy0,
Integral niter )

Applies a spectral collocation method based on Chebyshev nodes for solving differential equations.

This function utilizes a spectral collocation method based on Chebyshev nodes to solve a differential equation over an interval from x = x0 to x = x0 + h. The solution starts from the initial conditions y(x0) = y0 and ‘y’(x0) = dy0. In each iteration of the spectral collocation method, the number of Chebyshev nodes used increases, starting frominfo.niniand doubling until info.nmaxis reached or the desired tolerance is met. Theniterparameter keeps track of the number of iterations and is used to retrieve pre-computed differentiation matrices and Chebyshev nodes from theinfo` object.

Template Parameters
SolverInfoA type containing pre-computed information for the solver, like differentiation matrices and Chebyshev nodes.
ScalarThe scalar type of the independent variable
YScalarThe scalar type of the dependent variable
IntegralThe integral type of the number of iterations
AllocatorAn allocator type
Parameters
infoSolverInfo object - Contains pre-computed information for the solver, like differentiation matrices and Chebyshev nodes.
x0The starting value of the independent variable.
hStep size for the spectral method.
y0Initial value of the dependent variable at x0.
dy0Initial derivative of the dependent variable at x0.
niterCounter for the number of iterations of the spectral collocation step performed.
allocAn allocator for the Eigen objects.
Returns
A tuple containing:
  1. Eigen::Vector<std::complex<Scalar>, Eigen::Dynamic, 1> - Numerical estimate of the solution at the end of the step, at x0 + h.
  2. Eigen::Vector<std::complex<Scalar>, Eigen::Dynamic, 1> - Numerical estimate of the derivative of the solution at the end of the step, at x0 + h.
  3. Eigen::VectorXd (real) - Chebyshev nodes used for the current iteration of the spectral collocation method, scaled to lie in the interval [x0, x0 + h].

Definition at line 393 of file chebyshev.hpp.

References eval(), gamma(), omega(), scale(), and square().

◆ sqrt() [1/2]

template<typename T , require_not_floating_point< T > * = nullptr>
auto riccati::sqrt ( T && x)
inline

Definition at line 186 of file utils.hpp.

◆ sqrt() [2/2]

template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto riccati::sqrt ( T x)
inline

Definition at line 181 of file utils.hpp.

◆ square() [1/2]

template<typename T , require_not_floating_point< T > * = nullptr>
auto riccati::square ( T && x)
inline

Definition at line 196 of file utils.hpp.

◆ square() [2/2]

template<typename T , require_floating_point_or_complex< T > * = nullptr>
auto riccati::square ( T x)
inline

Definition at line 191 of file utils.hpp.

◆ step()

template<typename SolverInfo , typename Scalar , typename Vec >
auto riccati::step ( SolverInfo & info,
Scalar xi,
Scalar xf,
std::complex< Scalar > yi,
std::complex< Scalar > dyi,
Scalar eps,
Scalar epsilon_h,
Scalar init_stepsize,
Vec && x_eval,
bool hard_stop = false )
inline

Solves the differential equation y'' + 2gy' + w^2y = 0 over a given interval.

This function solves the differential equation on the interval (xi, xf), starting from the initial conditions y(xi) = yi and y'(xi) = dyi. It keeps the residual of the ODE below eps, and returns an interpolated solution (dense output) at the points specified in x_eval.

Template Parameters
SolverInfoType of the solver info object containing differentiation matrices, etc.
ScalarNumeric scalar type, typically float or double.
VecType of the vector for dense output values, should match Scalar type.
SolverInfoType of the solver info object containing differentiation matrices, etc.
ScalarNumeric scalar type, typically float or double.
VecType of the vector for dense output values.
Parameters
[in]infoSolverInfo object containing necessary information for the solver.
[in]xiStarting value of the independent variable.
[in]xfEnding value of the independent variable.
[in]yiInitial value of the dependent variable at xi.
[in]dyiInitial derivative of the dependent variable at xi.
[in]epsRelative tolerance for the local error of both Riccati and Chebyshev type steps.
[in]epsilon_hRelative tolerance for choosing the stepsize of Riccati steps.
[in]init_stepsizeinitial stepsize for the integration
[in]x_evalList of x-values where the solution is to be interpolated (dense output) and returned.
[in]hard_stopIf true, forces the solver to have a potentially smaller last stepsize to stop exactly at xf.
Returns
A tuple containing multiple elements representing the results of the ODE solving process: 0. std::vector<Scalar>: A vector containing the x-values at which the solution was evaluated or interpolated. These values correspond to the points in the interval [xi, xf] and include the points specified in x_eval if dense output was requested.
  1. std::vector<std::complex<Scalar>>: A vector of complex numbers representing the solution y(x) of the differential equation at each x-value from the corresponding vector of x-values.
  2. std::vector<std::complex<Scalar>>: A vector of complex numbers representing the derivative of the solution, y'(x), at each x-value from the corresponding vector of x-values.
  3. std::vector<int>: A vector indicating the success status of the solver at each step. Each element corresponds to a step in the solver process, where a value of 1 indicates success, and 0 indicates failure.
  4. std::vector<int>: A vector indicating the type of step taken at each point in the solution process. Each element corresponds to a step in the solver process, where a value of 1 indicates an oscillatory step and 0 indicates a non-oscillatory step.
  5. std::vector<Scalar>: A vector containing the phase angle at each step of the solution process if relevant. This is especially applicable for oscillatory solutions and may not be used for all types of differential equations.
  6. Eigen::Matrix<std::complex<Scalar>, -1, 1>: A vector containing the interpolated solution at the specified x_eval The function returns these vectors encapsulated in a standard tuple, providing comprehensive information about the solution process, including where the solution was evaluated, the values and derivatives of the solution at those points, success status of each step, type of each step, and phase angles where applicable.

Definition at line 273 of file step.hpp.

References riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::alloc_, riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::chebyshev_, choose_nonosc_stepsize(), choose_osc_stepsize(), compile_size_v, riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::Dn(), eval(), gamma(), get_slice(), interpolate(), matrix(), nonosc_step(), omega(), osc_step(), pi(), scale(), riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::xn(), and riccati::SolverInfo< OmegaFun, GammaFun, Scalar_, Integral_, Allocator >::xp().

◆ to_arena() [1/2]

template<typename T , typename Expr >
auto riccati::to_arena ( arena_allocator< T, arena_alloc > & arena,
const Expr & expr )
inlinenoexcept

Definition at line 155 of file arena_matrix.hpp.

◆ to_arena() [2/2]

template<typename Expr >
auto riccati::to_arena ( dummy_allocator & arena,
const Expr & expr )
inlinenoexcept

Definition at line 167 of file arena_matrix.hpp.

References eval().

◆ zero_like() [1/2]

template<typename T , require_not_floating_point< T > * = nullptr>
auto riccati::zero_like ( const T & x)
inline

Definition at line 226 of file utils.hpp.

◆ zero_like() [2/2]

template<typename T , require_floating_point_or_complex< T > * = nullptr>
T riccati::zero_like ( T x)
inlineconstexpr

Definition at line 221 of file utils.hpp.

Variable Documentation

◆ compile_size_v

template<typename T >
Eigen::Index riccati::compile_size_v = std::decay_t<T>::RowsAtCompileTime * std::decay_t<T>::ColsAtCompileTime
constexpr

Definition at line 13 of file utils.hpp.