Template Numerical Library version\ main:978fbd2
Loading...
Searching...
No Matches
TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true > Struct Template Reference
Inheritance diagram for TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >:
Collaboration diagram for TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >:

Public Types

using IndexType = std::size_t
 Type for indexing.
 
using RealType = GetValueType_t< Vector >
 Type of floating-point arithemtics.
 
using SolverMonitorType = SolverMonitor
 Type of object used for monitoring the convergence.
 
using ValueType = typename VectorType::ValueType
 
using VectorType = Vector
 
- Public Types inherited from TNL::Solvers::ODE::StaticExplicitSolver< GetValueType_t< Vector >, std::size_t >
using IndexType
 Indexing type.
 
using RealType
 Type of the floating-point arithmetics or static vector.
 

Public Member Functions

__cuda_callable__ ODESolver ()
 Default constructor.
 
__cuda_callable__ ODESolver (const ODESolver &solver)
 
__cuda_callable__ RealType getAdaptivity () const
 Getter of the parameter controlling the adaptive choice of the integration time step.
 
__cuda_callable__ Method & getMethod ()
 Gives reference to the underlying method.
 
__cuda_callable__ const Method & getMethod () const
 Gives constant reference to the underlying method.
 
void __cuda_callable__ init (const VectorType &u)
 Setup auxiliary vectors of the solver.
 
template<typename RHSFunction , typename... Params>
void __cuda_callable__ iterate (VectorType &u, RealType &time, RealType &tau, RHSFunction &&f, Params &&... params)
 Performs one iteration of the solver.
 
void __cuda_callable__ reset ()
 This method is just for consistency with the ODE solver for dynamic vectors.
 
__cuda_callable__ void setAdaptivity (const RealType &adaptivity)
 Setter of the parameter controlling the adaptive choice of the integration time step.
 
bool setup (const Config::ParameterContainer &parameters, const String &prefix="")
 Method for setup of the explicit solver based on configuration parameters.
 
template<typename RHSFunction , typename... Params>
__cuda_callable__ bool solve (VectorType &u, RHSFunction &&f, Params &&... params)
 Solve ODE given by a lambda function.
 
- Public Member Functions inherited from TNL::Solvers::ODE::StaticExplicitSolver< GetValueType_t< Vector >, std::size_t >
__cuda_callable__ StaticExplicitSolver ()=default
 Default constructor.
 
bool __cuda_callable__ checkNextIteration ()
 Checks if the solver is allowed to do the next iteration.
 
__cuda_callable__ const RealTypegetMaxTau () const
 Getter of maximal value of the time step.
 
__cuda_callable__ const RealTypegetStopTime () const
 Getter of the time where the evolution computation shall by stopped.
 
__cuda_callable__ const RealTypegetTau () const
 Getter of the time step used for the computation.
 
__cuda_callable__ const RealTypegetTime () const
 Getter of the current time of the evolution computed by the solver.
 
__cuda_callable__ void setMaxTau (const RealType &maxTau)
 Setter of maximal value of the time step.
 
__cuda_callable__ void setStopTime (const RealType &stopTime)
 Setter of the time where the evolution computation shall by stopped.
 
__cuda_callable__ void setTau (const RealType &tau)
 Setter of the time step used for the computation.
 
__cuda_callable__ void setTestingMode (bool testingMode)
 
__cuda_callable__ void setTime (const RealType &t)
 Settter of the current time of the evolution computed by the solver.
 
bool setup (const Config::ParameterContainer &parameters, const std::string &prefix="")
 Method for setup of the iterative solver based on configuration parameters.
 
- Public Member Functions inherited from TNL::Solvers::StaticIterativeSolver< GetValueType_t< Vector >, std::size_t >
__cuda_callable__ StaticIterativeSolver ()=default
 Default constructor.
 
__cuda_callable__ bool checkConvergence ()
 Checks whether the convergence occurred already.
 
__cuda_callable__ bool checkNextIteration ()
 Checks if the solver is allowed to do the next iteration.
 
__cuda_callable__ const GetValueType_t< Vector > & getConvergenceResidue () const
 Gets the the convergence threshold.
 
__cuda_callable__ const GetValueType_t< Vector > & getDivergenceResidue () const
 Gets the limit for the divergence criterion.
 
__cuda_callable__ const std::size_tgetIterations () const
 Gets the number of iterations performed by the solver so far.
 
__cuda_callable__ const std::size_tgetMaxIterations () const
 Gets the maximal number of iterations the solver is allowed to perform.
 
__cuda_callable__ const std::size_tgetMinIterations () const
 Gets the minimal number of iterations the solver is supposed to do.
 
__cuda_callable__ const GetValueType_t< Vector > & getResidue () const
 Gets the residue reached at the current iteration.
 
__cuda_callable__ bool nextIteration ()
 Proceeds to the next iteration.
 
__cuda_callable__ void resetIterations ()
 Sets the the number of the current iterations to zero.
 
__cuda_callable__ void setConvergenceResidue (const GetValueType_t< Vector > &convergenceResidue)
 Sets the threshold for the convergence.
 
__cuda_callable__ void setDivergenceResidue (const GetValueType_t< Vector > &divergenceResidue)
 Sets the residue limit for the divergence criterion.
 
__cuda_callable__ void setMaxIterations (const std::size_t &maxIterations)
 Sets the maximal number of iterations the solver is allowed to perform.
 
__cuda_callable__ void setMinIterations (const std::size_t &minIterations)
 Sets the minimal number of iterations the solver is supposed to do.
 
__cuda_callable__ void setResidue (const GetValueType_t< Vector > &residue)
 Sets the residue reached at the current iteration.
 
bool setup (const Config::ParameterContainer &parameters, const std::string &prefix="")
 Method for setup of the iterative solver based on configuration parameters.
 

Static Public Member Functions

static void configSetup (Config::ConfigDescription &config, const String &prefix="")
 Static method for setup of configuration parameters.
 
static constexpr bool isStatic ()
 
- Static Public Member Functions inherited from TNL::Solvers::ODE::StaticExplicitSolver< GetValueType_t< Vector >, std::size_t >
static void configSetup (Config::ConfigDescription &config, const std::string &prefix="")
 This method defines configuration entries for setup of the iterative solver.
 
- Static Public Member Functions inherited from TNL::Solvers::StaticIterativeSolver< GetValueType_t< Vector >, std::size_t >
static void configSetup (Config::ConfigDescription &config, const std::string &prefix="")
 This method defines configuration entries for setup of the iterative solver.
 

Static Public Attributes

static constexpr int Stages = Method::getStages()
 

Protected Attributes

RealType adaptivity = 0.00001
 
std::array< VectorType, Stages > k_vectors
 
VectorType kAux
 
Method method
 
- Protected Attributes inherited from TNL::Solvers::ODE::StaticExplicitSolver< GetValueType_t< Vector >, std::size_t >
RealType maxTau
 
bool stopOnSteadyState
 
RealType stopTime
 
RealType tau
 
bool testingMode
 
RealType time
 
- Protected Attributes inherited from TNL::Solvers::StaticIterativeSolver< GetValueType_t< Vector >, std::size_t >
GetValueType_t< Vector > convergenceResidue
 
std::size_t currentIteration
 
GetValueType_t< Vector > currentResidue
 
GetValueType_t< Vector > divergenceResidue
 
std::size_t maxIterations
 
std::size_t minIterations
 
std::size_t refreshRate
 

Member Typedef Documentation

◆ SolverMonitorType

template<typename Method , typename Vector , typename SolverMonitor >
using TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::SolverMonitorType = SolverMonitor

Type of object used for monitoring the convergence.

Can be TNL::Solvers::IterativeSolverMonitor.

Member Function Documentation

◆ configSetup()

template<typename Method , typename Vector , typename SolverMonitor >
static void TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::configSetup ( Config::ConfigDescription & config,
const String & prefix = "" )
static

Static method for setup of configuration parameters.

Parameters
configis the config description.
prefixis the prefix of the configuration parameters for this solver.

◆ getAdaptivity()

template<typename Method , typename Vector , typename SolverMonitor >
__cuda_callable__ RealType TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::getAdaptivity ( ) const
inline

Getter of the parameter controlling the adaptive choice of the integration time step.

Returns
the current value of the parameter controlling the adaptive choice of integration time step.

◆ getMethod() [1/2]

template<typename Method , typename Vector , typename SolverMonitor >
__cuda_callable__ Method & TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::getMethod ( )

Gives reference to the underlying method.

Returns
reference to the underlying method.

◆ getMethod() [2/2]

template<typename Method , typename Vector , typename SolverMonitor >
__cuda_callable__ const Method & TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::getMethod ( ) const

Gives constant reference to the underlying method.

Returns
constant reference to the underlying method.

◆ init()

template<typename Method , typename Vector , typename SolverMonitor >
void __cuda_callable__ TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::init ( const VectorType & u)

Setup auxiliary vectors of the solver.

This method is supposed to be called before the first call of the method iterate. It is not necessary to call this method before the method solve is used.

Parameters
uthis parameter is only for consistency with the ODE solver for dynamic vectors.

◆ iterate()

template<typename Method , typename Vector , typename SolverMonitor >
template<typename RHSFunction , typename... Params>
void __cuda_callable__ TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::iterate ( VectorType & u,
RealType & time,
RealType & tau,
RHSFunction && f,
Params &&... params )

Performs one iteration of the solver.

This method can be used for hybrid solvers which combine various ODE solvers. Otherwise, use of solve is recommended. Before the first call of this method, the method init has to be called.

Template Parameters
RHSFunctionis type of a lambda function representing the right-hand side of the ODE system.
Paramsare the parameters which are supposed to be passed to the lambda function f.
Parameters
uis a variable/static vector representing the solution of the ODE system at current time.
timeis the current time of the evolution. The variable is increased by tau.
tauis the current time step. It can be changed by the solver if the adaptive time step control is used.
fis the lambda function representing the right-hand side of the ODE system. The definition of the lambda function is the same as in the method solve.
paramsare the parameters which are supposed to be passed to the lambda function f.
Example
#include <iostream>
#include <TNL/Solvers/ODE/ODESolver.h>
#include <TNL/Solvers/ODE/Methods/Euler.h>
using Real = double;
int
main( int argc, char* argv[] )
{
const Real final_time = 10.0;
const Real output_time_step = 0.25;
Real next_output_time = TNL::min( output_time_step, final_time );
Real tau = 0.001;
Real time = 0.0;
Real current_tau = tau;
Vector u = 0.0;
ODESolver solver;
solver.init( u );
while( time < final_time ) {
auto f = []( const Real& t, const Real& current_tau, const Vector& u, Vector& fu )
{
fu = t * sin( t );
};
solver.iterate( u, time, current_tau, f );
if( time >= next_output_time ) {
std::cout << time << " " << u[ 0 ] << std::endl;
next_output_time += output_time_step;
if( next_output_time > final_time )
next_output_time = final_time;
}
if( time + current_tau > next_output_time )
current_tau = next_output_time - time;
else
current_tau = tau;
}
}
Vector with constant size.
Definition StaticVector.h:19
Definition Real.h:14
T endl(T... args)
constexpr ResultType min(const T1 &a, const T2 &b)
This function returns minimum of two numbers.
Definition Math.h:23
First order Euler method.
Definition Euler.h:15
Integrator or solver of systems of ordinary differential equations.
Definition ODESolver.h:61

◆ setAdaptivity()

template<typename Method , typename Vector , typename SolverMonitor >
__cuda_callable__ void TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::setAdaptivity ( const RealType & adaptivity)
inline

Setter of the parameter controlling the adaptive choice of the integration time step.

The smaller the parameter is the smaller the integration time step tends to be. Reasonable values for this parameters are approximately from interval \( [10^{-12},10^{-2}] \).

Parameters
adaptivitynew value of the parameter controlling the adaptive choice of integration time step.

◆ setup()

template<typename Method , typename Vector , typename SolverMonitor >
bool TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::setup ( const Config::ParameterContainer & parameters,
const String & prefix = "" )

Method for setup of the explicit solver based on configuration parameters.

Parameters
parametersis the container for configuration parameters.
prefixis the prefix of the configuration parameters for this solver.
Returns
true if the parameters where parsed successfully.
false if the method did not succeed to read the configuration parameters.

◆ solve()

template<typename Method , typename Vector , typename SolverMonitor >
template<typename RHSFunction , typename... Params>
__cuda_callable__ bool TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, true >::solve ( VectorType & u,
RHSFunction && f,
Params &&... params )

Solve ODE given by a lambda function.

Template Parameters
RHSFunctionis type of a lambda function representing the right-hand side of the ODE system. The definition of the lambda function reads as:
auto f = [=] ( const Real& t, const Real& tau, const VectorType& u, VectorType& fu ) {...}
where t is the current time of the evolution, tau is the current time step, u is the solution at the current time, fu is variable/static vector into which the lambda function is suppsed to evaluate the function \( f(t, \vec x) \) at the current time \( t \).
Parameters
uis a variable/static vector representing the solution of the ODE system at current time.
fis the lambda function representing the right-hand side of the ODE system.
paramsare the parameters which are supposed to be passed to the lambda function f. This is due to the fact that the CUDA compiler does not allow nested lambda functions: "An extended __host__ __device__ lambda cannot be defined inside an extended __host__ __device__ lambda expression".
Returns
true if steady state solution has been reached, false` otherwise.
Example
#include <iostream>
#include <TNL/Solvers/ODE/ODESolver.h>
#include <TNL/Solvers/ODE/Methods/Euler.h>
using Real = double;
int
main( int argc, char* argv[] )
{
const Real final_t = 10.0;
const Real tau = 0.001;
const Real output_time_step = 0.25;
ODESolver solver;
solver.setTau( tau );
solver.setTime( 0.0 );
Vector u = 0.0;
while( solver.getTime() < final_t ) {
solver.setStopTime( TNL::min( solver.getTime() + output_time_step, final_t ) );
auto f = []( const Real& t, const Real& tau, const Vector& u, Vector& fu )
{
fu = t * sin( t );
};
solver.solve( u, f );
std::cout << solver.getTime() << " " << u[ 0 ] << std::endl;
}
}
#include <iostream>
#include <fstream>
#include <TNL/Containers/Vector.h>
#include <TNL/Algorithms/parallelFor.h>
#include <TNL/Containers/StaticVector.h>
#include <TNL/Solvers/ODE/ODESolver.h>
#include <TNL/Solvers/ODE/Methods/Euler.h>
using Real = double;
template< typename Device >
void
solveParallelODEs( const char* file_name )
{
const Real final_t = 50.0;
const Real tau = 0.001;
const Real output_time_step = 0.005;
const Real sigma_min( 10.0 ), rho_min( 15.0 ), beta_min( 1.0 );
const int parametric_steps = 4;
const Real sigma_step = 30.0 / ( parametric_steps - 1 );
const Real rho_step = 21.0 / ( parametric_steps - 1 );
const Real beta_step = 15.0 / ( parametric_steps - 1 );
const int output_time_steps = ceil( final_t / output_time_step ) + 1;
const int results_size( output_time_steps * parametric_steps * parametric_steps * parametric_steps );
auto results_view = results.getView();
auto f = [ = ] __cuda_callable__( const Real& t,
const Real& tau,
const StaticVector& u,
StaticVector& fu,
const Real& sigma_i,
const Real& rho_j,
const Real& beta_k )
{
const Real& x = u[ 0 ];
const Real& y = u[ 1 ];
const Real& z = u[ 2 ];
fu[ 0 ] = sigma_i * ( y - x );
fu[ 1 ] = rho_j * x - y - x * z;
fu[ 2 ] = -beta_k * z + x * y;
};
auto solve = [ = ] __cuda_callable__( const MultiIndex& idx ) mutable
{
const Real sigma_i = sigma_min + idx[ 0 ] * sigma_step;
const Real rho_j = rho_min + idx[ 1 ] * rho_step;
const Real beta_k = beta_min + idx[ 2 ] * beta_step;
ODESolver solver;
solver.setTau( tau );
solver.setTime( 0.0 );
StaticVector u( 1.0, 1.0, 1.0 );
int time_step( 1 );
results_view[ ( idx[ 0 ] * parametric_steps + idx[ 1 ] ) * parametric_steps + idx[ 2 ] ] = u;
while( time_step < output_time_steps ) {
solver.setStopTime( TNL::min( solver.getTime() + output_time_step, final_t ) );
solver.solve( u, f, sigma_i, rho_j, beta_k );
const int k =
( ( time_step++ * parametric_steps + idx[ 0 ] ) * parametric_steps + idx[ 1 ] ) * parametric_steps + idx[ 2 ];
results_view[ k ] = u;
}
};
const MultiIndex begin = { 0, 0, 0 };
const MultiIndex end = { parametric_steps, parametric_steps, parametric_steps };
file.open( file_name, std::ios::out );
for( int sigma_idx = 0; sigma_idx < parametric_steps; sigma_idx++ )
for( int rho_idx = 0; rho_idx < parametric_steps; rho_idx++ )
for( int beta_idx = 0; beta_idx < parametric_steps; beta_idx++ ) {
Real sigma = sigma_min + sigma_idx * sigma_step;
Real rho = rho_min + rho_idx * rho_step;
Real beta = beta_min + beta_idx * beta_step;
file << "# sigma " << sigma << " rho " << rho << " beta " << beta << std::endl;
for( int i = 0; i < output_time_steps - 1; i++ ) {
int offset = ( ( i * parametric_steps + sigma_idx ) * parametric_steps + rho_idx ) * parametric_steps + beta_idx;
auto u = results.getElement( offset );
file << u[ 0 ] << " " << u[ 1 ] << " " << u[ 2 ] << std::endl;
}
file << std::endl;
}
}
int
main( int argc, char* argv[] )
{
if( argc != 2 ) {
std::cout << "Usage: " << argv[ 0 ] << " <path to output directory>" << std::endl;
return EXIT_FAILURE;
}
TNL::String file_name( argv[ 1 ] );
file_name += "/StaticODESolver-LorenzParallelExample-result.out";
solveParallelODEs< TNL::Devices::Host >( file_name.getString() );
#ifdef __CUDACC__
solveParallelODEs< TNL::Devices::Cuda >( file_name.getString() );
#endif
}
#define __cuda_callable__
Definition Macros.h:49
Array with constant size.
Definition StaticArray.h:20
Vector extends Array with algebraic operations.
Definition Vector.h:36
Class for managing strings.
Definition String.h:30
std::enable_if_t< std::is_integral_v< Begin > &&std::is_integral_v< End > > parallelFor(const Begin &begin, const End &end, typename Device::LaunchConfiguration launch_config, Function f, FunctionArgs... args)
Parallel for-loop function for 1D range specified with integral values.
Definition parallelFor.h:41
T open(T... args)

.


The documentation for this struct was generated from the following file: