allpix::RungeKutta
Class to perform arbitrary Runge-Kutta integration. More…
#include <runge_kutta.h>
Public Classes
Name | |
---|---|
class | Step Utility type to return both the value and the error at every step. |
Public Types
Name | |
---|---|
using std::function< Eigen::Matrix< T, D, 1 >(T, Eigen::Matrix< T, D, 1 >)> | StepFunction Stepping function to integrate a single step of the equations. |
Public Functions
Name | |
---|---|
RungeKutta(Eigen::Matrix< T, S+2, S > tableau, StepFunction function, T step_size, Eigen::Matrix< T, D, 1 > initial_y, T initial_t =0) Construct a Runge-Kutta integrator. |
|
void | setTimeStep(T step_size) Changes the time step. |
T | getTimeStep() const Return the time step. |
void | setValue(Eigen::Matrix< T, D, 1 > y) Changes the current value during integration. |
Eigen::Matrix< T, D, 1 > | getValue() const Get the value to integrate. |
Eigen::Matrix< T, D, 1 > | getError() const Get the total integration error. |
T | getTime() const Get the time during integration. |
void | advanceTime(double t) Advance the time of the integration. |
Step | step() Execute a single time step of the integration. |
Step | step(int amount) Execute multiple time steps of the integration. |
Detailed Description
template <typename T ,
int S,
int D =3>
class allpix::RungeKutta;
Class to perform arbitrary Runge-Kutta integration.
Class can be provided a Runge-Kutta tableau (optionally with an error function), together with the dimension of the equations and a step function to integrate a step of the equation. Both the result, error and timestep can be retrieved and changed during the integration.
Public Types Documentation
using StepFunction
using allpix::RungeKutta< T, S, D >::StepFunction = std::function<Eigen::Matrix<T, D, 1>(T, Eigen::Matrix<T, D, 1>)>;
Stepping function to integrate a single step of the equations.
Public Functions Documentation
function RungeKutta
inline RungeKutta(
Eigen::Matrix< T, S+2, S > tableau,
StepFunction function,
T step_size,
Eigen::Matrix< T, D, 1 > initial_y,
T initial_t =0
)
Construct a Runge-Kutta integrator.
Parameters:
- tableau One of the possible Runge-Kutta tables (see allpix::tableau should be preferred)
- function Step function to perform integration
- step_size Time step of the integration
- initial_y Start values of the vector to perform integration on
- initial_t Initial time at the start of the integration
function setTimeStep
inline void setTimeStep(
T step_size
)
Changes the time step.
Parameters:
- step_size New time step of the integration
function getTimeStep
inline T getTimeStep() const
Return the time step.
Return: Current time step of the integration
function setValue
inline void setValue(
Eigen::Matrix< T, D, 1 > y
)
Changes the current value during integration.
Note: Can be used to add additional processes during the integration
function getValue
inline Eigen::Matrix< T, D, 1 > getValue() const
Get the value to integrate.
Return: Current value
function getError
inline Eigen::Matrix< T, D, 1 > getError() const
Get the total integration error.
Return: Total integrated error
function getTime
inline T getTime() const
Get the time during integration.
Return: Current time
function advanceTime
inline void advanceTime(
double t
)
Advance the time of the integration.
Parameters:
- t Time step to advance the integration by
function step
inline Step step()
Execute a single time step of the integration.
Return: Combination of the current value and the error in this single step
function step
inline Step step(
int amount
)
Execute multiple time steps of the integration.
Parameters:
- amount Number of steps to combine
Return: Combination of the current value and the total error in all the steps
Updated on 2024-12-13 at 08:31:36 +0000