Lorenz equations field.
It is a FOODIE integrand class concrete extension.
The Lorenz' equations system [1] is a non linear system of pure ODEs that retains a reasonable-complex behaviour: such a system, for a certain parameters-region exhibits a chaotic dynamics useful for testing FOODIE solvers.
The Lorenz' ODEs system can be written as:
$$\begin{matrix} U_t = R(U) \\ U = \begin{bmatrix} v_1 \\ v_2 \\ v_3 \end{bmatrix}\;\;\; R(U) = \begin{bmatrix} \sigma (v_2-v_1) \\ v_1(\rho - v_3) -v_2 \\ v_1 v_2 - \beta v_3 \end{bmatrix} \end{matrix}$$
The parameters set is constant and it is here selected as:
$$\begin{matrix} \sigma = 10 \\ \rho = 28 \\ \beta = \frac{8}{3} \end{matrix}$$
These values are chaos-inducing thus they magnify the eventual numerical inaccuracies of FOODIE solvers, see [2].
[1] Deterministic Nonperiodic Flow, Lorenz E.N., Journal of the Atmospheric Sciences, 1963, vol. 20, pp. 130–141, doi: http://dx.doi.org/10.1175/1520-0469(1963)020<0130:DNF>2.0.CO;2
[2] Scientific software design: the object-oriented way, Rouson, Damian, Jim Xia, and Xiaofeng Xu, Cambridge University Press, 2011
State variables are organized as an array (rank 1) of reals of dims elements, in this case 3 elements.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=R_P), | private, | dimension(:), allocatable | :: | U | Integrand (state) variables, [1:dims]. |
||
real(kind=R_P), | private | :: | beta | = | 0._R_P | Lorenz \(\beta\). |
|
integer(kind=I_P), | private | :: | dims | = | 0 | Space dimensions. |
|
real(kind=R_P), | private, | dimension(:,:), allocatable | :: | previous | Previous time steps states, [1:dims,1:steps]. |
||
real(kind=R_P), | private | :: | rho | = | 0._R_P | Lorenz \(\rho\). |
|
real(kind=R_P), | private | :: | sigma | = | 0._R_P | Lorenz \(\sigma\). |
|
integer(kind=I_P), | private | :: | steps | = | 0 | Number of time steps stored. |
Lorenz + Lorenz operator.
Add two Lorenz fields.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(in) | :: | lhs | Left hand side. |
||
class(integrand), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Lorenz = Lorenz.
Assign one Lorenz field to another.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(inout) | :: | lhs | Left hand side. |
||
class(integrand), | intent(in) | :: | rhs | Right hand side. |
Lorenz = real.
Assign one real to a Lorenz field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(inout) | :: | lhs | Left hand side. |
||
real(kind=R_P), | intent(in) | :: | rhs | Right hand side. |
Init field.
Construct an initialized Lorenz field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(inout) | :: | self | Lorenz field. |
||
real(kind=R_P), | intent(in), | dimension(:) | :: | initial_state | Initial state of Lorenz field vector. |
|
real(kind=R_P), | intent(in) | :: | sigma | Lorenz \(\sigma\). |
||
real(kind=R_P), | intent(in) | :: | rho | Lorenz \(\rho\). |
||
real(kind=R_P), | intent(in) | :: | beta | Lorenz \(\beta\). |
||
integer(kind=I_P), | intent(in), | optional | :: | steps | Time steps stored. |
Lorenz * Lorenz operator.
Multiply a lorenz field by another one.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(in) | :: | lhs | Left hand side. |
||
class(integrand), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Lorenz * real operator.
Multiply a Lorenz field by a real scalar.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(in) | :: | lhs | Left hand side. |
||
real(kind=R_P), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Local error.
Estimate local truncation error between 2 lorenz approximations.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(in) | :: | lhs | Left hand side. |
||
class(integrand), | intent(in) | :: | rhs | Right hand side. |
Error estimation.
Extract Lorenz field.
Get a previous time step.
Extract previous time solution of Lorenz field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(in) | :: | self | Lorenz field. |
||
integer(kind=I_P), | intent(in) | :: | n | Time level. |
Previous time solution of Lorenz field.
Real * Lorenz operator.
Multiply a real scalar by a Lorenz field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=R_P), | intent(in) | :: | lhs | Left hand side. |
||
class(lorenz), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Lorenz - Lorenz.
Subtract two Lorenz fields.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(in) | :: | lhs | Left hand side. |
||
class(integrand), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Time derivative, residuals function.
Time derivative of Lorenz field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(in) | :: | self | Lorenz field. |
||
real(kind=R_P), | intent(in), | optional | :: | t | Time. |
Lorenz field time derivative.
Update previous time steps.
Update previous time steps.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(lorenz), | intent(inout) | :: | self | Lorenz field. |
||
class(integrand), | intent(in), | optional | :: | filter | Filter field displacement. |
|
real(kind=R_P), | intent(in), | optional | :: | weights(:) | Weights for filtering the steps. |
type, extends(integrand) :: lorenz
!< Lorenz equations field.
!<
!< It is a FOODIE integrand class concrete extension.
!<
!<### Lorenz ODEs system
!<The Lorenz' equations system [1] is a non linear system of pure ODEs that retains a reasonable-complex behaviour: such a
!<system, for a certain parameters-region exhibits a chaotic dynamics useful for testing FOODIE solvers.
!<
!<The Lorenz' ODEs system can be written as:
!<
!<$$\begin{matrix}
!< U_t = R(U) \\
!< U = \begin{bmatrix}
!< v_1 \\
!< v_2 \\
!< v_3
!< \end{bmatrix}\;\;\;
!< R(U) = \begin{bmatrix}
!< \sigma (v_2-v_1) \\
!< v_1(\rho - v_3) -v_2 \\
!< v_1 v_2 - \beta v_3
!< \end{bmatrix}
!<\end{matrix}$$
!<
!<The parameters set is constant and it is here selected as:
!<
!<$$\begin{matrix}
!< \sigma = 10 \\
!< \rho = 28 \\
!< \beta = \frac{8}{3}
!<\end{matrix}$$
!<
!<These values are chaos-inducing thus they magnify the eventual numerical inaccuracies of FOODIE solvers, see [2].
!<
!<#### Bibliography
!<
!<[1] *Deterministic Nonperiodic Flow*, Lorenz E.N., Journal of the Atmospheric Sciences, 1963, vol. 20, pp. 130--141,
!<doi: http://dx.doi.org/10.1175/1520-0469(1963)020<0130:DNF>2.0.CO;2
!<
!<[2] *Scientific software design: the object-oriented way*, Rouson, Damian, Jim Xia, and Xiaofeng Xu,
!<Cambridge University Press, 2011
!<
!<#### State variables organization
!< State variables are organized as an array (rank 1) of reals of *dims* elements, in this case 3 elements.
private
integer(I_P) :: dims=0 !< Space dimensions.
integer(I_P) :: steps=0 !< Number of time steps stored.
real(R_P), dimension(:), allocatable :: U !< Integrand (state) variables, [1:dims].
real(R_P), dimension(:,:), allocatable :: previous !< Previous time steps states, [1:dims,1:steps].
real(R_P) :: sigma=0._R_P !< Lorenz \(\sigma\).
real(R_P) :: rho=0._R_P !< Lorenz \(\rho\).
real(R_P) :: beta=0._R_P !< Lorenz \(\beta\).
contains
! auxiliary methods
procedure, pass(self), public :: init !< Init field.
procedure, pass(self), public :: output !< Extract Lorenz field.
! ADT integrand deferred methods
procedure, pass(self), public :: t => dLorenz_dt !< Time derivative, residuals function.
procedure, pass(lhs), public :: local_error => lorenz_local_error !< Local error.
procedure, pass(self), public :: update_previous_steps !< Update previous time steps.
procedure, pass(self), public :: previous_step !< Get a previous time step.
procedure, pass(lhs), public :: integrand_multiply_integrand => lorenz_multiply_lorenz !< Lorenz * Lorenz operator.
procedure, pass(lhs), public :: integrand_multiply_real => lorenz_multiply_real !< Lorenz * real operator.
procedure, pass(rhs), public :: real_multiply_integrand => real_multiply_lorenz !< Real * Lorenz operator.
procedure, pass(lhs), public :: add => add_lorenz !< Lorenz + Lorenz operator.
procedure, pass(lhs), public :: sub => sub_lorenz !< Lorenz - Lorenz.
procedure, pass(lhs), public :: assign_integrand => lorenz_assign_lorenz !< Lorenz = Lorenz.
procedure, pass(lhs), public :: assign_real => lorenz_assign_real !< Lorenz = real.
endtype lorenz