template<typename Method, typename Vector, typename
SolverMonitor = IterativeSolverMonitor< typename Vector::RealType, typename Vector::IndexType >, bool IsStatic = IsStaticArrayType< Vector >()>
struct TNL::Solvers::ODE::ODESolver< Method, Vector, SolverMonitor, IsStatic >
Integrator or solver of systems of ordinary differential equations.
\( \vec u( 0 ) = \vec u_{ini} \).
The vector \( \vec u(t) \) can be represented using different types of containers, depending on the size and nature of the ODE system:
The method, which is supposed to be used by the solver, is represented by the template parameter Method.
3#include <TNL/Containers/Vector.h>
4#include <TNL/Algorithms/parallelFor.h>
5#include <TNL/Containers/StaticVector.h>
6#include <TNL/Solvers/ODE/ODESolver.h>
7#include <TNL/Solvers/ODE/Methods/Euler.h>
14template<
typename Device >
16solveParallelODEs(
const char* file_name )
25 const Real final_t = 50.0;
26 const Real tau = 0.001;
27 const Real output_time_step = 0.005;
31 const Real sigma_min( 10.0 ), rho_min( 15.0 ), beta_min( 1.0 );
32 const int parametric_steps = 4;
33 const Real sigma_step = 30.0 / ( parametric_steps - 1 );
34 const Real rho_step = 21.0 / ( parametric_steps - 1 );
35 const Real beta_step = 15.0 / ( parametric_steps - 1 );
36 const int output_time_steps = ceil( final_t / output_time_step ) + 1;
40 const int results_size( output_time_steps * parametric_steps * parametric_steps * parametric_steps );
42 auto results_view = results.getView();
48 const StaticVector& u,
54 const Real& x = u[ 0 ];
55 const Real& y = u[ 1 ];
56 const Real& z = u[ 2 ];
57 fu[ 0 ] = sigma_i * ( y - x );
58 fu[ 1 ] = rho_j * x - y - x * z;
59 fu[ 2 ] = -beta_k * z + x * y;
67 const Real sigma_i = sigma_min + idx[ 0 ] * sigma_step;
68 const Real rho_j = rho_min + idx[ 1 ] * rho_step;
69 const Real beta_k = beta_min + idx[ 2 ] * beta_step;
75 solver.setTime( 0.0 );
76 StaticVector u( 1.0, 1.0, 1.0 );
78 results_view[ ( idx[ 0 ] * parametric_steps + idx[ 1 ] ) * parametric_steps + idx[ 2 ] ] = u;
82 while( time_step < output_time_steps ) {
83 solver.setStopTime( TNL::min( solver.getTime() + output_time_step, final_t ) );
85 solver.solve( u, f, sigma_i, rho_j, beta_k );
88 ( ( time_step++ * parametric_steps + idx[ 0 ] ) * parametric_steps + idx[ 1 ] ) * parametric_steps + idx[ 2 ];
89 results_view[ k ] = u;
96 const MultiIndex begin = { 0, 0, 0 };
97 const MultiIndex end = { parametric_steps, parametric_steps, parametric_steps };
103 file.
open( file_name, std::ios::out );
104 for(
int sigma_idx = 0; sigma_idx < parametric_steps; sigma_idx++ )
105 for(
int rho_idx = 0; rho_idx < parametric_steps; rho_idx++ )
106 for(
int beta_idx = 0; beta_idx < parametric_steps; beta_idx++ ) {
107 Real sigma = sigma_min + sigma_idx * sigma_step;
108 Real rho = rho_min + rho_idx * rho_step;
109 Real beta = beta_min + beta_idx * beta_step;
110 file <<
"# sigma " << sigma <<
" rho " << rho <<
" beta " << beta <<
std::endl;
111 for(
int i = 0; i < output_time_steps - 1; i++ ) {
112 int offset = ( ( i * parametric_steps + sigma_idx ) * parametric_steps + rho_idx ) * parametric_steps + beta_idx;
113 auto u = results.getElement( offset );
114 file << u[ 0 ] <<
" " << u[ 1 ] <<
" " << u[ 2 ] <<
std::endl;
122main(
int argc,
char* argv[] )
129 file_name +=
"/StaticODESolver-LorenzParallelExample-result.out";
131 solveParallelODEs< TNL::Devices::Host >( file_name.getString() );
133 solveParallelODEs< TNL::Devices::Cuda >( file_name.getString() );
#define __cuda_callable__
Definition Macros.h:49
Array with constant size.
Definition StaticArray.h:20
Vector with constant size.
Definition StaticVector.h:19
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
First order Euler method.
Definition Euler.h:15
Integrator or solver of systems of ordinary differential equations.
Definition ODESolver.h:61
2#include <TNL/FileName.h>
3#include <TNL/Containers/Vector.h>
4#include <TNL/Solvers/ODE/ODESolver.h>
5#include <TNL/Solvers/ODE/Methods/Euler.h>
11template<
typename Device >
13solveHeatEquation(
const char* file_name )
17 using VectorView =
typename Vector::ViewType;
26 const Real final_t = 0.05;
27 const Real output_time_step = 0.005;
29 const Real h = 1.0 / ( n - 1 );
30 const Real tau = 0.1 * h * h;
31 const Real h_sqr_inv = 1.0 / ( h * h );
43 if( x >= 0.4 && x <= 0.6 )
49 file.
open( file_name, std::ios::out );
50 write( file, u, n, h, (
Real) 0.0 );
59 solver.setTime( 0.0 );
66 while( solver.getTime() < final_t ) {
67 solver.setStopTime( TNL::min( solver.getTime() + output_time_step, final_t ) );
69 auto f = [ = ]
__cuda_callable__( Index i,
const VectorView& u, VectorView& fu )
mutable
71 if( i == 0 || i == n - 1 )
74 fu[ i ] = h_sqr_inv * ( u[ i - 1 ] - 2.0 * u[ i ] + u[ i + 1 ] );
78 auto time_stepping = [ = ](
const Real& t,
const Real& tau,
const VectorView& u, VectorView& fu )
85 solver.solve( u, time_stepping );
86 write( file, u, n, h, solver.getTime() );
92main(
int argc,
char* argv[] )
99 file_name +=
"/ODESolver-HeatEquationExample-result.out";
101 solveHeatEquation< TNL::Devices::Host >( file_name.getString() );
103 solveHeatEquation< TNL::Devices::Cuda >( file_name.getString() );