A Modular Arbitrary-Order Ocean-Atmosphere Model -- Stochastic implementation
Functions/Subroutines | Variables
rk2_mtv_integrator Module Reference

Module with the MTV rk2 integration routines. More...

Functions/Subroutines

subroutine, public init_integrator
 Subroutine to initialize the MTV rk2 integrator. More...
 
subroutine init_noise
 Routine to initialize the noise vectors and buffers. More...
 
subroutine init_g
 Routine to initialize the G term. More...
 
subroutine compg (y)
 Routine to actualize the G term based on the state y of the MTV system. More...
 
subroutine, public step (y, t, dt, dtn, res, tend)
 Routine to perform an integration step (Heun algorithm) of the MTV system. The incremented time is returned. More...
 
subroutine, public full_step (y, t, dt, dtn, res)
 Routine to perform an integration step (Heun algorithm) of the full stochastic system. The incremented time is returned. More...
 

Variables

real(kind=8), dimension(:), allocatable buf_y1
 
real(kind=8), dimension(:), allocatable buf_f0
 
real(kind=8), dimension(:), allocatable buf_f1
 Integration buffers. More...
 
real(kind=8), dimension(:), allocatable dw
 
real(kind=8), dimension(:), allocatable dwmult
 Standard gaussian noise buffers. More...
 
real(kind=8), dimension(:), allocatable dwar
 
real(kind=8), dimension(:), allocatable dwau
 
real(kind=8), dimension(:), allocatable dwor
 
real(kind=8), dimension(:), allocatable dwou
 Standard gaussian noise buffers. More...
 
real(kind=8), dimension(:), allocatable anoise
 
real(kind=8), dimension(:), allocatable noise
 Additive noise term. More...
 
real(kind=8), dimension(:), allocatable noisemult
 Multiplicative noise term. More...
 
real(kind=8), dimension(:), allocatable g
 G term of the MTV tendencies. More...
 
real(kind=8), dimension(:), allocatable buf_g
 Buffer for the G term computation. More...
 
logical mult
 Logical indicating if the sigma1 matrix must be computed for every state change. More...
 
logical q1fill
 Logical indicating if the matrix Q1 is non-zero. More...
 
logical compute_mult
 Logical indicating if the Gaussian noise for the multiplicative noise must be computed. More...
 
real(kind=8), parameter sq2 = sqrt(2.D0)
 Hard coded square root of 2. More...
 

Detailed Description

Module with the MTV rk2 integration routines.

Remarks
This module actually contains the Heun algorithm routines.

Function/Subroutine Documentation

subroutine rk2_mtv_integrator::compg ( real(kind=8), dimension(0:ndim), intent(in)  y)
private

Routine to actualize the G term based on the state y of the MTV system.

Parameters
yState of the MTV system

Definition at line 105 of file rk2_MTV_integrator.f90.

105  REAL(KIND=8), DIMENSION(0:ndim), INTENT(IN) :: y
106 
107  g=htot
108  CALL sparse_mul2_k(ltot,y,buf_g)
109  g=g+buf_g
110  CALL sparse_mul3(btot,y,y,buf_g)
111  g=g+buf_g
112  CALL sparse_mul4(mtot,y,y,y,buf_g)
113  g=g+buf_g
subroutine, public rk2_mtv_integrator::full_step ( real(kind=8), dimension(0:ndim), intent(in)  y,
real(kind=8), intent(inout)  t,
real(kind=8), intent(in)  dt,
real(kind=8), intent(in)  dtn,
real(kind=8), dimension(0:ndim), intent(out)  res 
)

Routine to perform an integration step (Heun algorithm) of the full stochastic system. The incremented time is returned.

Parameters
yInitial point.
tActual integration time
dtIntegration timestep.
dtnStochastoc integration timestep (normally square-root of dt).
resFinal point after the step.

Definition at line 170 of file rk2_MTV_integrator.f90.

170  REAL(KIND=8), DIMENSION(0:ndim), INTENT(IN) :: y
171  REAL(KIND=8), INTENT(INOUT) :: t
172  REAL(KIND=8), INTENT(IN) :: dt,dtn
173  REAL(KIND=8), DIMENSION(0:ndim), INTENT(OUT) :: res
174  CALL stoch_atm_res_vec(dwar)
175  CALL stoch_atm_unres_vec(dwau)
176  CALL stoch_oc_res_vec(dwor)
177  CALL stoch_oc_unres_vec(dwou)
178  anoise=(q_ar*dwar+q_au*dwau+q_or*dwor+q_ou*dwou)*dtn
179  CALL sparse_mul3(aotensor,y,y,buf_f0)
180  buf_y1 = y+dt*buf_f0+anoise
181  CALL sparse_mul3(aotensor,buf_y1,buf_y1,buf_f1)
182  res=y+0.5*(buf_f0+buf_f1)*dt+anoise
183  t=t+dt
type(coolist), dimension(:), allocatable, public aotensor
- Tensor representation of the tendencies.
subroutine rk2_mtv_integrator::init_g ( )
private

Routine to initialize the G term.

Definition at line 97 of file rk2_MTV_integrator.f90.

97  INTEGER :: allocstat
98  ALLOCATE(g(0:ndim), buf_g(0:ndim), stat=allocstat)
99  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
integer ndim
Number of variables (dimension of the model)
Definition: params.f90:85
Statistics accumulators.
Definition: stat.f90:14
subroutine, public rk2_mtv_integrator::init_integrator ( )

Subroutine to initialize the MTV rk2 integrator.

Definition at line 50 of file rk2_MTV_integrator.f90.

50  INTEGER :: allocstat
51 
52  CALL init_ss_integrator ! Initialize the uncoupled resolved dynamics
53 
54  ALLOCATE(buf_y1(0:ndim),buf_f0(0:ndim),buf_f1(0:ndim),stat=allocstat)
55  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
56 
57  buf_y1=0.d0
58  buf_f1=0.d0
59  buf_f0=0.d0
60 
61  print*, 'Initializing the integrator ...'
62  CALL init_sigma(mult,q1fill)
63  CALL init_noise
64  CALL init_g
integer ndim
Number of variables (dimension of the model)
Definition: params.f90:85
Statistics accumulators.
Definition: stat.f90:14
subroutine rk2_mtv_integrator::init_noise ( )
private

Routine to initialize the noise vectors and buffers.

Definition at line 69 of file rk2_MTV_integrator.f90.

69  INTEGER :: allocstat
70  ALLOCATE(dw(0:ndim), dwmult(0:ndim), stat=allocstat)
71  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
72 
73  ALLOCATE(dwar(0:ndim),dwau(0:ndim),dwor(0:ndim),dwou(0:ndim),stat=allocstat)
74  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
75 
76  ALLOCATE(anoise(0:ndim), noise(0:ndim), noisemult(0:ndim), stat=allocstat)
77  IF (allocstat /= 0) stop "*** Not enough memory ! ***"
78 
79  dw=0.d0
80  dwmult=0.d0
81 
82  dwar=0.d0
83  dwor=0.d0
84  dwau=0.d0
85  dwou=0.d0
86 
87  anoise=0.d0
88  noise=0.d0
89  noisemult=0.d0
90 
91  compute_mult=((q1fill).OR.(mult))
92 
integer ndim
Number of variables (dimension of the model)
Definition: params.f90:85
Statistics accumulators.
Definition: stat.f90:14
subroutine, public rk2_mtv_integrator::step ( real(kind=8), dimension(0:ndim), intent(in)  y,
real(kind=8), intent(inout)  t,
real(kind=8), intent(in)  dt,
real(kind=8), intent(in)  dtn,
real(kind=8), dimension(0:ndim), intent(out)  res,
real(kind=8), dimension(0:ndim), intent(out)  tend 
)

Routine to perform an integration step (Heun algorithm) of the MTV system. The incremented time is returned.

Parameters
yInitial point.
tActual integration time
dtIntegration timestep.
dtnStochastic integration timestep (normally square-root of dt).
resFinal point after the step.
tendPartial or full tendencies used to perform the step (used for debugging).

Definition at line 124 of file rk2_MTV_integrator.f90.

124  REAL(KIND=8), DIMENSION(0:ndim), INTENT(IN) :: y
125  REAL(KIND=8), INTENT(INOUT) :: t
126  REAL(KIND=8), INTENT(IN) :: dt,dtn
127  REAL(KIND=8), DIMENSION(0:ndim), INTENT(OUT) :: res,tend
128 
129  CALL compg(y)
130 
131  CALL stoch_atm_res_vec(dwar)
132  CALL stoch_oc_res_vec(dwor)
133  anoise=q_ar*dwar+q_or*dwor
134  CALL stoch_vec(dw)
135  IF (compute_mult) CALL stoch_vec(dwmult)
136  noise(1:ndim)=matmul(sig2,dw(1:ndim))
137  IF ((mult).and.(mod(t,mnuti)<dt)) CALL compute_mult_sigma(y)
138  IF (compute_mult) noisemult(1:ndim)=matmul(sig1,dwmult(1:ndim))
139 
140  CALL tendencies(t,y,buf_f0)
141  buf_y1 = y+dt*(buf_f0+g)+(anoise+sq2*(noise+noisemult))*dtn
142 
143  buf_f1=g
144  CALL compg(buf_y1)
145  g=0.5*(g+buf_f1)
146 
147  IF ((mult).and.(mod(t,mnuti)<dt)) CALL compute_mult_sigma(buf_y1)
148  IF (compute_mult) THEN
149  buf_f1(1:ndim)=matmul(sig1,dwmult(1:ndim))
150  noisemult(1:ndim)=0.5*(noisemult(1:ndim)+buf_f1(1:ndim))
151  ENDIF
152 
153 
154  CALL tendencies(t,buf_y1,buf_f1)
155  buf_f0=0.5*(buf_f0+buf_f1)
156  res=y+dt*(buf_f0+g)+(anoise+sq2*(noise+noisemult))*dtn
157  ! tend=G+sq2*(noise+noisemult)/dtn
158  tend=sq2*noisemult/dtn
159  t=t+dt
160 
integer ndim
Number of variables (dimension of the model)
Definition: params.f90:85

Variable Documentation

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::anoise
private

Definition at line 33 of file rk2_MTV_integrator.f90.

33  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: anoise,noise !< Additive noise term
real(kind=8), dimension(:), allocatable rk2_mtv_integrator::buf_f0
private

Definition at line 30 of file rk2_MTV_integrator.f90.

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::buf_f1
private

Integration buffers.

Definition at line 30 of file rk2_MTV_integrator.f90.

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::buf_g
private

Buffer for the G term computation.

Definition at line 36 of file rk2_MTV_integrator.f90.

36  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: buf_g !< Buffer for the G term computation
real(kind=8), dimension(:), allocatable rk2_mtv_integrator::buf_y1
private

Definition at line 30 of file rk2_MTV_integrator.f90.

30  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: buf_y1,buf_f0,buf_f1 !< Integration buffers
logical rk2_mtv_integrator::compute_mult
private

Logical indicating if the Gaussian noise for the multiplicative noise must be computed.

Definition at line 40 of file rk2_MTV_integrator.f90.

40  LOGICAL :: compute_mult !< Logical indicating if the Gaussian noise for the multiplicative noise must be computed
real(kind=8), dimension(:), allocatable rk2_mtv_integrator::dw
private

Definition at line 31 of file rk2_MTV_integrator.f90.

31  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: dw, dwmult !< Standard gaussian noise buffers
real(kind=8), dimension(:), allocatable rk2_mtv_integrator::dwar
private

Definition at line 32 of file rk2_MTV_integrator.f90.

32  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: dwar,dwau,dwor,dwou !< Standard gaussian noise buffers
real(kind=8), dimension(:), allocatable rk2_mtv_integrator::dwau
private

Definition at line 32 of file rk2_MTV_integrator.f90.

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::dwmult
private

Standard gaussian noise buffers.

Definition at line 31 of file rk2_MTV_integrator.f90.

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::dwor
private

Definition at line 32 of file rk2_MTV_integrator.f90.

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::dwou
private

Standard gaussian noise buffers.

Definition at line 32 of file rk2_MTV_integrator.f90.

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::g
private

G term of the MTV tendencies.

Definition at line 35 of file rk2_MTV_integrator.f90.

35  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: g !< G term of the MTV tendencies
logical rk2_mtv_integrator::mult
private

Logical indicating if the sigma1 matrix must be computed for every state change.

Definition at line 38 of file rk2_MTV_integrator.f90.

38  LOGICAL :: mult !< Logical indicating if the sigma1 matrix must be computed for every state change
real(kind=8), dimension(:), allocatable rk2_mtv_integrator::noise
private

Additive noise term.

Definition at line 33 of file rk2_MTV_integrator.f90.

real(kind=8), dimension(:), allocatable rk2_mtv_integrator::noisemult
private

Multiplicative noise term.

Definition at line 34 of file rk2_MTV_integrator.f90.

34  REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: noisemult !< Multiplicative noise term
logical rk2_mtv_integrator::q1fill
private

Logical indicating if the matrix Q1 is non-zero.

Definition at line 39 of file rk2_MTV_integrator.f90.

39  LOGICAL :: q1fill !< Logical indicating if the matrix Q1 is non-zero
real(kind=8), parameter rk2_mtv_integrator::sq2 = sqrt(2.D0)
private

Hard coded square root of 2.

Definition at line 42 of file rk2_MTV_integrator.f90.

42  REAL(KIND=8), PARAMETER :: sq2 = sqrt(2.d0) !< Hard coded square root of 2