Class WKBSolver
Defined in File wkbsolver.hpp
Inheritance Relationships
Derived Types
public WKBSolver1(Class WKBSolver1)public WKBSolver2(Class WKBSolver2)public WKBSolver3(Class WKBSolver3)
Class Documentation
-
class WKBSolver
Class to carry out WKB steps of varying orders.
Subclassed by WKBSolver1, WKBSolver2, WKBSolver3
Public Functions
-
WKBSolver()
-
Eigen::Matrix<std::complex<double>, 3, 2> step(std::complex<double> x0, std::complex<double> dx0, double t0, double h0, const Eigen::Matrix<std::complex<double>, 6, 1> &ws, const Eigen::Matrix<std::complex<double>, 6, 1> &gs, const Eigen::Matrix<std::complex<double>, 5, 1> &ws5, const Eigen::Matrix<std::complex<double>, 5, 1> &gs5)
Computes a WKB step of a given order and returns the solution and its local error estimate.
- Parameters:
x0[in] – value of the solution \(x(t)\) at the start of the step
dx0[in] – value of the derivative of the solution \(\frac{dx}{dt}\) at the start of the step
t0[in] – value of independent variable (time) at the start of the step
h0[in] – size of the step (solution will be given at t+t0)
ws[in] – vector of 6 evaulations of the frequency term at the Gauss-Lobatto nodes, this is necessary for Gauss-Lobatto integration, which in turn is needed to calculate the WKB series
gs[in] – vector of 6 evaluations of the friction term
ws5[in] – vector of 5 evaluations of the friction term at the nodes of 5th order Gauss-Lobatto quadrature, this is needed to compute the error on Gauss-Lobatto quadrature, and in turn the WKB series
gs5[in] – vector of 5 evaluations of the friction term
- Returns:
a matrix, whose rows are:
\( x, \dot{x}\) at t0+h0 as an nth order WKB estimate
\( \Delta_{\mathrm{trunc}}x, \Delta_{\mathrm{trunc}}\dot{x}\) at t0+h0, defined as the difference between an nth and (n-1)th order WKB estimate
\( \Delta_{\mathrm{int}}x, \Delta_{\mathrm{int}}\dot{x}\) at t0+h0, defined as the local error coming from those terms in the WKB series that involved numerical integrals
-
void dense_step(double t0, const std::list<double> &dots, std::list<std::complex<double>> &doxs, std::list<std::complex<double>> &dodxs)
Computes dense output at a set of timepoints within a step.
- Parameters:
t0[in] – value of independent variable (time) at the start of the step
dots[in] – sequence of timepoints at which dense output is to be generated
doxs[in, out] – dense output for the solution \(x(t)\)
dodxs[in, out] – dense output for the derivative of the solution \(\dot{x}\)
-
Eigen::Matrix<double, 6, 1> dense_weights_6(double t)
-
Eigen::Matrix<double, 6, 1> dense_weights_derivs_6(double t)
-
std::complex<double> dense_integrate(const Eigen::Matrix<double, 6, 1> &denseweights6, const Eigen::Matrix<std::complex<double>, 6, 1> &integrand6)
-
std::complex<double> dense_interpolate(const Eigen::Matrix<double, 6, 1> &denseweights6, const Eigen::Matrix<std::complex<double>, 6, 1> &integrand6)
Public Members
-
Eigen::Matrix<std::complex<double>, 7, 1> x_vdm
Protected Functions
-
void d1w1()
-
void d1w2()
-
void d1w3()
-
void d1w4()
-
void d1w5()
-
void d1w6()
-
void d2w1()
-
void d2w2()
-
void d2w3()
-
void d2w4()
-
void d2w5()
-
void d2w6()
-
void d3w1()
-
void d3w2()
-
void d3w3()
-
void d3w4()
-
void d3w5()
-
void d3w6()
-
void d4w1()
-
void d1g1()
-
void d1g2()
-
void d1g3()
-
void d1g4()
-
void d1g5()
-
void d1g6()
-
void d2g1()
-
void d2g2()
-
void d2g3()
-
void d2g4()
-
void d2g5()
-
void d2g6()
-
void d3g1()
-
void d1w2_5()
-
void d1w3_5()
-
void d1w4_5()
-
virtual void dds()
-
virtual void dsi()
-
virtual void dsf()
-
virtual void s()
-
void fp()
-
void fm()
-
void dfpi()
-
void dfmi()
-
void dfpf()
-
void dfmf()
-
void ddfp()
-
void ddfm()
-
void ap()
-
void am()
-
void bp()
-
void bm()
-
Eigen::Matrix<std::complex<double>, 2, 1> integrate(const Eigen::Matrix<std::complex<double>, 6, 1> &integrand6, const Eigen::Matrix<std::complex<double>, 5, 1> &integrand5)
Protected Attributes
-
Eigen::Matrix<double, 6, 1> glws6
-
Eigen::Matrix<double, 5, 1> glws5
-
Eigen::Matrix<double, 7, 1> d4w1_w
-
Eigen::Matrix<double, 6, 1> d1w1_w
-
Eigen::Matrix<double, 6, 1> d1w2_w
-
Eigen::Matrix<double, 6, 1> d1w3_w
-
Eigen::Matrix<double, 6, 1> d1w4_w
-
Eigen::Matrix<double, 6, 1> d1w5_w
-
Eigen::Matrix<double, 6, 1> d1w6_w
-
Eigen::Matrix<double, 6, 1> d2w1_w
-
Eigen::Matrix<double, 6, 1> d2w2_w
-
Eigen::Matrix<double, 6, 1> d2w3_w
-
Eigen::Matrix<double, 6, 1> d2w4_w
-
Eigen::Matrix<double, 6, 1> d2w5_w
-
Eigen::Matrix<double, 6, 1> d2w6_w
-
Eigen::Matrix<double, 6, 1> d3w1_w
-
Eigen::Matrix<double, 6, 1> d3w2_w
-
Eigen::Matrix<double, 6, 1> d3w3_w
-
Eigen::Matrix<double, 6, 1> d3w4_w
-
Eigen::Matrix<double, 6, 1> d3w5_w
-
Eigen::Matrix<double, 6, 1> d3w6_w
-
Eigen::Matrix<double, 6, 1> d1g1_w
-
Eigen::Matrix<double, 6, 1> d1g6_w
-
Eigen::Matrix<double, 6, 1> d2g1_w
-
Eigen::Matrix<double, 6, 1> d2g6_w
-
Eigen::Matrix<double, 6, 1> d3g1_w
-
Eigen::Matrix<double, 5, 1> d1w2_5_w
-
Eigen::Matrix<double, 5, 1> d1w3_5_w
-
Eigen::Matrix<double, 5, 1> d1w4_5_w
-
Eigen::Matrix<std::complex<double>, 7, 1> ws7_
-
Eigen::Matrix<std::complex<double>, 6, 1> ws_
-
Eigen::Matrix<std::complex<double>, 6, 1> gs_
-
Eigen::Matrix<std::complex<double>, 5, 1> ws5_
-
Eigen::Matrix<std::complex<double>, 5, 1> gs5_
-
std::complex<double> d1w1_
-
std::complex<double> d1w2_
-
std::complex<double> d1w3_
-
std::complex<double> d1w4_
-
std::complex<double> d1w5_
-
std::complex<double> d1w6_
-
std::complex<double> d2w1_
-
std::complex<double> d2w2_
-
std::complex<double> d2w3_
-
std::complex<double> d2w4_
-
std::complex<double> d2w5_
-
std::complex<double> d2w6_
-
std::complex<double> d3w1_
-
std::complex<double> d3w2_
-
std::complex<double> d3w3_
-
std::complex<double> d3w4_
-
std::complex<double> d3w5_
-
std::complex<double> d3w6_
-
std::complex<double> d4w1_
-
std::complex<double> d1g1_
-
std::complex<double> d1g2_
-
std::complex<double> d1g3_
-
std::complex<double> d1g4_
-
std::complex<double> d1g5_
-
std::complex<double> d1g6_
-
std::complex<double> d2g1_
-
std::complex<double> d2g2_
-
std::complex<double> d2g3_
-
std::complex<double> d2g4_
-
std::complex<double> d2g5_
-
std::complex<double> d2g6_
-
std::complex<double> d3g1_
-
std::complex<double> d1w2_5_
-
std::complex<double> d1w3_5_
-
std::complex<double> d1w4_5_
-
Eigen::Matrix<std::complex<double>, 6, 1> dws_
-
Eigen::Matrix<std::complex<double>, 6, 1> dgs_
-
Eigen::Matrix<std::complex<double>, 6, 1> d2ws_
-
Eigen::Matrix<std::complex<double>, 6, 1> d2gs_
-
Eigen::Matrix<std::complex<double>, 6, 1> d3ws_
-
Eigen::Matrix<std::complex<double>, 5, 1> dws5_
-
Eigen::Matrix<std::complex<double>, 1, 4> dds_
-
Eigen::Matrix<std::complex<double>, 1, 4> dsi_
-
Eigen::Matrix<std::complex<double>, 1, 4> dsf_
-
Eigen::Matrix<std::complex<double>, 1, 4> s_
-
Eigen::Matrix<std::complex<double>, 1, 4> s_error
-
std::complex<double> fp_
-
std::complex<double> fm_
-
std::complex<double> dfpi_
-
std::complex<double> dfmi_
-
std::complex<double> dfpf_
-
std::complex<double> dfmf_
-
std::complex<double> ddfp_
-
std::complex<double> ddfm_
-
std::complex<double> ap_
-
std::complex<double> am_
-
std::complex<double> bp_
-
std::complex<double> bm_
-
double h
-
std::complex<double> x
-
std::complex<double> dx
-
std::complex<double> ddx
-
int order_
-
std::complex<double> err_fp
-
std::complex<double> err_fm
-
std::complex<double> err_dfp
-
std::complex<double> err_dfm
-
std::list<std::complex<double>> doxs
-
std::list<std::complex<double>> dodxs
-
std::list<std::complex<double>> dows
-
Eigen::Matrix<std::complex<double>, 1, 4> dense_s_
-
Eigen::Matrix<std::complex<double>, 1, 4> dense_ds_
-
Eigen::Matrix<std::complex<double>, 1, 4> dense_ds_i
-
std::complex<double> dense_ap_
-
std::complex<double> dense_am_
-
std::complex<double> dense_bp_
-
std::complex<double> dense_bm_
-
Eigen::Matrix<double, 6, 6> integ_vandermonde
-
Eigen::Matrix<double, 6, 6> interp_vandermonde
-
WKBSolver()