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

Public Types

using DeviceType = typename Vector::DeviceType
 Device where the solver is supposed to be executed.
 
using DofVectorType = TNL::Containers::Vector< RealType, DeviceType, IndexType >
 Alias for type of unknown variable \( \vec x \).
 
using IndexType = typename Vector::IndexType
 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 Vector::ValueType
 
using VectorType = Vector
 Type of unknown variable \( \vec x \).
 
- Public Types inherited from TNL::Solvers::ODE::ExplicitSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
using IndexType
 Indexing type.
 
using RealType
 Type of the floating-point arithmetics.
 
using SolverMonitorType
 Type of the monitor of the convergence of the solver.
 
- Public Types inherited from TNL::Solvers::IterativeSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
using SolverMonitorType
 Type of an object used for monitoring of the convergence.
 

Public Member Functions

 ODESolver ()
 Default constructor.
 
__cuda_callable__ RealType getAdaptivity () const
 Getter of the parameter controlling the adaptive choice of the integration time step.
 
Method & getMethod ()
 Gives reference to the underlying method.
 
const Method & getMethod () const
 Gives constant reference to the underlying method.
 
void init (const VectorType &u)
 Setup auxiliary vectors of the solver.
 
template<typename RHSFunction , typename... Params>
void iterate (VectorType &u, RealType &time, RealType &tau, RHSFunction &&f, Params &&... params)
 Performs one iteration of the solver.
 
void reset ()
 Resets the solver.
 
__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>
bool solve (VectorType &u, RHSFunction &&f, Params &&... params)
 Solve ODE given by a lambda function.
 
- Public Member Functions inherited from TNL::Solvers::ODE::ExplicitSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
 ExplicitSolver ()=default
 Default constructor.
 
bool checkNextIteration ()
 Checks if the solver is allowed to do the next iteration.
 
const RealTypegetMaxTau () const
 Getter of maximal value of the time step.
 
const RealTypegetStopTime () const
 Getter of the time where the evolution computation shall by stopped.
 
const RealTypegetTau () const
 Getter of the time step used for the computation.
 
const RealTypegetTime () const
 Getter of the current time of the evolution computed by the solver.
 
void refreshSolverMonitor (bool force=false)
 This method refreshes the solver monitor.
 
void setMaxTau (const RealType &maxTau)
 Setter of maximal value of the time step.
 
void setStopTime (const RealType &stopTime)
 Setter of the time where the evolution computation shall by stopped.
 
void setTau (const RealType &tau)
 Setter of the time step used for the computation.
 
void setTestingMode (bool testingMode)
 
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::IterativeSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
 IterativeSolver ()=default
 Default constructor.
 
bool checkConvergence ()
 Checks whether the convergence occurred already.
 
bool checkNextIteration ()
 Checks if the solver is allowed to do the next iteration.
 
const Vector::RealType & getConvergenceResidue () const
 Gets the the convergence threshold.
 
const Vector::RealType & getDivergenceResidue () const
 Gets the limit for the divergence criterion.
 
const Vector::IndexType & getIterations () const
 Gets the number of iterations performed by the solver so far.
 
const Vector::IndexType & getMaxIterations () const
 Gets the maximal number of iterations the solver is allowed to perform.
 
const Vector::IndexType & getMinIterations () const
 Gets the minimal number of iterations the solver is supposed to do.
 
const Vector::RealType & getResidue () const
 Gets the residue reached at the current iteration.
 
bool nextIteration ()
 Proceeds to the next iteration.
 
void resetIterations ()
 Sets the the number of the current iterations to zero.
 
void setConvergenceResidue (const Vector::RealType &convergenceResidue)
 Sets the threshold for the convergence.
 
void setDivergenceResidue (const Vector::RealType &divergenceResidue)
 Sets the residue limit for the divergence criterion.
 
void setMaxIterations (const Vector::IndexType &maxIterations)
 Sets the maximal number of iterations the solver is allowed to perform.
 
void setMinIterations (const Vector::IndexType &minIterations)
 Sets the minimal number of iterations the solver is supposed to do.
 
void setRefreshRate (const Vector::IndexType &refreshRate)
 Sets the refresh rate (in milliseconds) for the solver monitor.
 
void setResidue (const Vector::RealType &residue)
 Sets the residue reached at the current iteration.
 
void setSolverMonitor (SolverMonitorType &solverMonitor)
 Sets the solver monitor object.
 
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::ExplicitSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
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::IterativeSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
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::ExplicitSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
RealType maxTau
 
bool stopOnSteadyState
 
RealType stopTime
 
RealType tau
 
bool testingMode
 
RealType time
 
- Protected Attributes inherited from TNL::Solvers::IterativeSolver< Vector::RealType, Vector::IndexType, SolverMonitor >
Vector::RealType convergenceResidue
 
Vector::IndexType currentIteration
 
Vector::RealType currentResidue
 
Vector::RealType divergenceResidue
 
Vector::IndexType maxIterations
 
Vector::IndexType minIterations
 
Vector::IndexType refreshRate
 
std::ofstream residualHistoryFile
 
std::string residualHistoryFileName
 
SolverMonitorsolverMonitor
 

Member Typedef Documentation

◆ DofVectorType

template<typename Method , typename Vector , typename SolverMonitor >
using TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::DofVectorType = TNL::Containers::Vector< RealType, DeviceType, IndexType >

Alias for type of unknown variable \( \vec x \).

Note, VectorType can be TNL::Containers::VectorView but DofVectorType is always TNL::Containers::Vector.

◆ SolverMonitorType

template<typename Method , typename Vector , typename SolverMonitor >
using TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::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, false >::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, false >::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 >
Method & TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::getMethod ( )

Gives reference to the underlying method.

Returns
reference to the underlying method.

◆ getMethod() [2/2]

template<typename Method , typename Vector , typename SolverMonitor >
const Method & TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::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 TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::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. Also this methods neeeds to be called everytime the size of u changes.

Parameters
uis a variable/dynamic vector representing the solution of the ODE system at current time.

◆ iterate()

template<typename Method , typename Vector , typename SolverMonitor >
template<typename RHSFunction , typename... Params>
void TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::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

◆ reset()

template<typename Method , typename Vector , typename SolverMonitor >
void TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::reset ( )

Resets the solver.

This method frees memory allocated by the solver. If it is called, the method init has to be called before the next call of the method iterate.

◆ setAdaptivity()

template<typename Method , typename Vector , typename SolverMonitor >
__cuda_callable__ void TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::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, false >::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>
bool TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, false >::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 ) {...}
Vector VectorType
Type of unknown variable .
Definition ODESolver.h:276
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/FileName.h>
#include <TNL/Containers/Vector.h>
#include <TNL/Solvers/ODE/ODESolver.h>
#include <TNL/Solvers/ODE/Methods/Euler.h>
#include <TNL/Solvers/IterativeSolverMonitor.h>
#include "write.h"
using Real = double;
using Index = int;
template< typename Device >
void
solveHeatEquation( const char* file_name )
{
using VectorView = typename Vector::ViewType;
/***
* Parameters of the discertisation
*/
const Real final_t = 0.05;
const Real output_time_step = 0.005;
const Index n = 41;
const Real h = 1.0 / ( n - 1 );
const Real tau = 0.001 * h * h;
const Real h_sqr_inv = 1.0 / ( h * h );
/***
* Initial condition
*/
Vector u( n );
u.forAllElements(
[ = ] __cuda_callable__( Index i, Real & value )
{
const Real x = i * h;
if( x >= 0.4 && x <= 0.6 )
value = 1.0;
else
value = 0.0;
} );
file.open( file_name, std::ios::out );
write( file, u, n, h, (Real) 0.0 );
/***
* Setup monitor for the ODE solver.
*/
using IterativeSolverMonitorType = TNL::Solvers::IterativeSolverMonitor< Real, Index >;
IterativeSolverMonitorType monitor;
TNL::Solvers::SolverMonitorThread monitorThread( monitor );
monitor.setRefreshRate( 10 ); // refresh rate in miliseconds
monitor.setVerbose( 1 );
monitor.setStage( "ODE solver stage:" );
/***
* Setup of the solver
*/
ODESolver solver;
solver.setTau( tau );
solver.setTime( 0.0 );
solver.setSolverMonitor( monitor );
/***
* Time loop
*/
while( solver.getTime() < final_t ) {
solver.setStopTime( TNL::min( solver.getTime() + output_time_step, final_t ) );
auto f = [ = ] __cuda_callable__( Index i, const VectorView& u, VectorView& fu ) mutable
{
if( i == 0 || i == n - 1 ) // boundary nodes -> boundary conditions
fu[ i ] = 0.0;
else // interior nodes -> approximation of the second derivative
fu[ i ] = h_sqr_inv * ( u[ i - 1 ] - 2.0 * u[ i ] + u[ i + 1 ] );
};
auto time_stepping = [ = ]( const Real& t, const Real& tau, const VectorView& u, VectorView& fu )
{
};
solver.solve( u, time_stepping );
write( file, u, n, h, solver.getTime() ); // write the current state to a file
}
monitor.stopMainLoop();
}
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 += "/ODESolver-HeatEquationExampleWithMonitor-result.out";
solveHeatEquation< TNL::Devices::Host >( file_name.getString() );
#ifdef __CUDACC__
solveHeatEquation< TNL::Devices::Cuda >( file_name.getString() );
#endif
}
#define __cuda_callable__
Definition Macros.h:49
Vector extends Array with algebraic operations.
Definition Vector.h:36
A RAII wrapper for launching the SolverMonitor's main loop in a separate thread.
Definition SolverMonitor.h:133
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: