Files
CoolProp/include/ODEIntegrators.h
Ian Bell 636fba7093 fix(docs): resolve all 162 doxygen warnings (#2708)
- Remove 13 obsolete Doxyfile tags (TCL_SUBST, CLANG_ASSISTED_PARSING,
  CLANG_OPTIONS, COLS_IN_ALPHA_INDEX, HTML_TIMESTAMP, FORMULA_TRANSPARENT,
  LATEX_SOURCE_CODE, PERL_PATH, CLASS_DIAGRAMS, MSCGEN_PATH, DOT_FONTNAME,
  DOT_FONTSIZE, DOT_TRANSPARENT)
- Fix unclosed \f$ inline math blocks in TransportRoutines.h and VLERoutines.h
- Fix spurious \f] inside display math block in TransportRoutines.h
- Fix @param name mismatches (matrix→mat, mole_fractions→mf, x0→x, etc.)
- Add missing @param docs to FlashRoutines, MixtureDerivatives,
  PhaseEnvelopeRoutines, ODEIntegrators, BicubicBackend, CoolProp.h
- Remove @return on void function in MixtureParameters.h
- Remove unresolvable \ref tags in ReducingFunctions.h and CoolPropLib.h
- Fix UNIFAC.h @param on no-arg function
- Remove ~25 duplicate /// @param blocks from PolyMath.cpp definitions
- Remove duplicate @param from IncompressibleBackend.cpp
- Clean up commented-out @param block leaking into calc_rhomass

Co-authored-by: Claude Sonnet 4.6 <noreply@anthropic.com>
2026-04-04 10:34:15 -04:00

41 lines
1.8 KiB
C++

#ifndef ODEINTEGRATORS_H
#define ODEINTEGRATORS_H
#include <vector>
namespace ODEIntegrators {
/// The abstract class defining the interface for the integrator routines
class AbstractODEIntegrator
{
public:
virtual std::vector<double> get_initial_array() const = 0;
virtual void pre_step_callback() = 0;
virtual void post_deriv_callback() = 0;
virtual void post_step_callback(double t, double h, std::vector<double>& x) = 0;
virtual bool premature_termination() = 0;
virtual void derivs(double t, std::vector<double>& x, std::vector<double>& f) = 0;
};
/**
@brief Use the adaptive Runge-Kutta integrator to integrate a system of differential equations
@param ode The ODE integrator instance that implements the required callbacks
@param tmin Starting value of the independent variable. ``t`` is in the closed range [``tmin``, ``tmax``]
@param tmax Ending value for the independent variable. ``t`` is in the closed range [``tmin``, ``tmax``]
@param hmin Minimum step size, something like 1e-5 usually is good. Don't make this too big or you may not be able to get a stable solution
@param hmax Maximum step size
@param eps_allowed Maximum absolute error of any CV per step allowed. Don't make this parameter too big or you may not be able to get a stable solution. Also don't make it too small because then you are going to run into truncation error.
@param step_relax The relaxation factor that is used in the step resizing algorithm. Should be less than 1.0; you can play with this parameter to improve the adaptive resizing, but should not be necessary.
*/
bool AdaptiveRK54(AbstractODEIntegrator& ode, double tmin, double tmax, double hmin, double hmax, double eps_allowed, double step_relax);
} // namespace ODEIntegrators
#endif