template<int Size_, class
Real>
class TNL::Solvers::ODE::StaticMerson< Containers::StaticVector< Size_, Real > >
Solver of ODEs with the first order of accuracy.
This solver is based on the Runge-Kutta-Merson method with adaptive choice of the integration time step for solving of ordinary differential equations having the following form:
\( \frac{d u}{dt} = f( t, u) \text{ on } (0,T) \)
\( u( 0 ) = u_{ini} \). It is supposed to be used when the unknown \( \vec x \in R^n \) is expressed by a Containers::StaticVector.
For problems where \( x\) is represented by floating-point number, see TNL::Solvers::ODE::StaticMerson. For problems where \( \vec x\) is represented by TNL::Containers::Vector or TNL::Containers::VectorView, see TNL::Solvers::ODE::Merson.
The following example demonstrates the use the solvers:
2#include <TNL/Containers/StaticVector.h>
3#include <TNL/Solvers/ODE/StaticEuler.h>
7int main(
int argc,
char* argv[] )
11 const Real final_t = 25.0;
12 const Real tau = 0.001;
13 const Real output_time_step = 0.01;
14 const Real sigma = 10.0;
15 const Real rho = 28.0;
16 const Real beta = 8.0 / 3.0;
20 solver.setTime( 0.0 );
21 Vector u( 1.0, 2.0, 3.0 );
22 while( solver.getTime() < final_t )
24 solver.setStopTime( TNL::min( solver.getTime() + output_time_step, final_t ) );
25 auto f = [=] (
const Real& t,
const Real& tau,
const Vector& u, Vector& fu ) {
26 const Real& x = u[ 0 ];
27 const Real& y = u[ 1 ];
28 const Real& z = u[ 2 ];
29 fu[ 0 ] = sigma * (y - x );
30 fu[ 1 ] = rho * x - y - x * z;
31 fu[ 2 ] = -beta * z + x * y;
Solver of ODEs with the first order of accuracy.
Definition StaticEuler.h:55
Since this variant of the Euler solver is static, it can be used even inside of GPU kernels and so combined with TNL::Algorithms::parallelFor as demonstrated by the following example:
3#include <TNL/Solvers/ODE/StaticEuler.h>
4#include <TNL/Containers/Vector.h>
5#include <TNL/Algorithms/parallelFor.h>
6#include <TNL/Containers/StaticArray.h>
11template<
typename Device >
12void solveParallelODEs(
const char* file_name )
16 const Real final_t = 50.0;
17 const Real tau = 0.001;
18 const Real output_time_step = 0.005;
19 const Real sigma_min( 10.0 ), rho_min( 15.0 ), beta_min( 1.0 );
20 const int parametric_steps = 4;
21 const Real sigma_step = 30.0 / ( parametric_steps - 1 );
22 const Real rho_step = 21.0 / ( parametric_steps - 1 );
23 const Real beta_step = 15.0 / ( parametric_steps - 1 );
24 const int output_time_steps = ceil( final_t / output_time_step ) + 1;
26 const int results_size( output_time_steps * parametric_steps * parametric_steps * parametric_steps );
28 auto results_view = results.getView();
30 const Real& sigma_i,
const Real& rho_j,
const Real& beta_k ) {
31 const Real& x = u[ 0 ];
32 const Real& y = u[ 1 ];
33 const Real& z = u[ 2 ];
34 fu[ 0 ] = sigma_i * (y - x );
35 fu[ 1 ] = rho_j * x - y - x * z;
36 fu[ 2 ] = -beta_k * z + x * y;
39 const Real sigma_i = sigma_min + i[ 0 ] * sigma_step;
40 const Real rho_j = rho_min + i[ 1 ] * rho_step;
41 const Real beta_k = beta_min + i[ 2 ] * beta_step;
45 solver.setTime( 0.0 );
46 Vector u( 1.0, 1.0, 1.0 );
48 results_view[ ( i[ 0 ] * parametric_steps + i[ 1 ] ) * parametric_steps + i[ 2 ] ] = u;
49 while( time_step < output_time_steps )
51 solver.setStopTime( TNL::min( solver.getTime() + output_time_step, final_t ) );
52 solver.solve( u, f, sigma_i, rho_j, beta_k );
53 const int idx = ( ( time_step++ * parametric_steps + i[ 0 ] ) * parametric_steps + i[ 1 ] ) * parametric_steps + i[ 2 ];
54 results_view[ idx ] = u;
57 const MultiIndex begin = { 0, 0, 0 };
58 const MultiIndex end = { parametric_steps, parametric_steps, parametric_steps };
59 TNL::Algorithms::parallelFor< Device >( begin, end,
solve );
62 file.
open( file_name, std::ios::out );
63 for(
int sigma_idx = 0; sigma_idx < parametric_steps; sigma_idx++ )
64 for(
int rho_idx = 0; rho_idx < parametric_steps; rho_idx++ )
65 for(
int beta_idx = 0; beta_idx < parametric_steps; beta_idx++ )
67 Real sigma = sigma_min + sigma_idx * sigma_step;
68 Real rho = rho_min + rho_idx * rho_step;
69 Real beta = beta_min + beta_idx * beta_step;
70 file <<
"# sigma " << sigma <<
" rho " << rho <<
" beta " << beta <<
std::endl;
71 for(
int i = 0; i < output_time_steps - 1; i++ )
73 int offset = ( ( i * parametric_steps + sigma_idx ) * parametric_steps + rho_idx ) * parametric_steps + beta_idx;
74 auto u = results.getElement( offset );
75 file << u[ 0 ] <<
" " << u[ 1 ] <<
" " << u[ 2 ] <<
std::endl;
81int main(
int argc,
char* argv[] )
84 file_name +=
"/StaticODESolver-LorenzParallelExample-result.out";
86 solveParallelODEs< TNL::Devices::Host >( file_name.getString() );
88 solveParallelODEs< TNL::Devices::Cuda >( file_name.getString() );
#define __cuda_callable__
Definition Macros.h:52
Array with constant size.
Definition StaticArray.h:23
Vector extends Array with algebraic operations.
Definition Vector.h:39
__cuda_callable__ bool solve(VectorType &u, RHSFunction &&f, Args... args)
Solve ODE given by a lambda function.
Class for managing strings.
Definition String.h:33
- Template Parameters
-
Size | is size of the ODE system. |
Real | is floating point number type, it is type of \( x \) in this case. |