accelerInt
v0.1
|
Classes | |
class | rhs_eval |
A wrapper class to evaluate the rhs function y' = f(y) stores the state variable, and provides to dydt. More... | |
Typedefs | |
typedef std::vector< double > | state_type |
state vector More... | |
typedef runge_kutta_fehlberg78< state_type, double > | stepper |
solver type More... | |
typedef controlled_runge_kutta< stepper > | controller |
controller type More... | |
Functions | |
void | initialize_solver (int) |
Initializes the solver. More... | |
void | cleanup_solver (int num_threads) |
Cleans up the created solvers. More... | |
const char * | solver_name () |
Returns a descriptive solver name. More... | |
void | init_solver_log () |
Initializes solver specific items for logging. More... | |
void | solver_log () |
Executes solver specific logging tasks. More... | |
void | intDriver (const int NUM, const double t, const double t_end, const double *pr_global, double *y_global) |
Integration driver for the CPU integrators. More... | |
Variables | |
std::vector< state_type * > | state_vectors |
State vector containers for boost. More... | |
std::vector< rhs_eval * > | evaluators |
RHS wrappers for boost. More... | |
std::vector< stepper * > | steppers |
Addaptive timesteppers. More... | |
std::vector< controller > | controllers |
ODE controllers. More... | |
typedef controlled_runge_kutta< stepper > rk78::controller |
controller type
Definition at line 35 of file rk78_typedefs.hpp.
typedef std::vector< double > rk78::state_type |
state vector
Definition at line 29 of file rk78_typedefs.hpp.
typedef runge_kutta_fehlberg78< state_type , double > rk78::stepper |
solver type
Definition at line 32 of file rk78_typedefs.hpp.
void rk78::cleanup_solver | ( | int | num_threads | ) |
Cleans up the created solvers.
num_threads | The number of OpenMP threads used |
Frees and cleans up allocated RK78 memory.
Definition at line 58 of file rk78_init.cpp.
void rk78::init_solver_log | ( | ) |
Initializes solver specific items for logging.
Initializes stepsize logging for stiffness measurement
Definition at line 86 of file rk78_init.cpp.
void rk78::initialize_solver | ( | int | num_threads | ) |
Initializes the solver.
num_threads | The number of OpenMP threads to use |
Definition at line 41 of file rk78_init.cpp.
void rk78::intDriver | ( | const int | NUM, |
const double | t, | ||
const double | t_end, | ||
const double * | pr_global, | ||
double * | y_global | ||
) |
Integration driver for the CPU integrators.
[in] | NUM | the number of IVPs to solve |
[in] | t | the current IVP time |
[in] | t_end | the time to integrate the IVP to |
[in] | pr_global | the pressure value for the IVPs |
[in,out] | y_global | the state vectors |
The integration driver for the RK78 solver
Definition at line 57 of file solver_rk78.cpp.
void rk78::solver_log | ( | ) |
Executes solver specific logging tasks.
Definition at line 96 of file rk78_init.cpp.
const char * rk78::solver_name | ( | ) |
Returns a descriptive solver name.
Definition at line 75 of file rk78_init.cpp.
std::vector< controller > rk78::controllers |
ODE controllers.
Definition at line 24 of file rk78_init.cpp.
std::vector< rhs_eval * > rk78::evaluators |
RHS wrappers for boost.
Definition at line 20 of file rk78_init.cpp.
std::vector< state_type * > rk78::state_vectors |
State vector containers for boost.
Definition at line 18 of file rk78_init.cpp.
std::vector<stepper*> rk78::steppers |
Addaptive timesteppers.
Definition at line 22 of file rk78_init.cpp.