New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
stprk3.F90 in NEMO/branches/2021/dev_r14393_HPC-03_Mele_Comm_Cleanup/src/SWE – NEMO

source: NEMO/branches/2021/dev_r14393_HPC-03_Mele_Comm_Cleanup/src/SWE/stprk3.F90 @ 14730

Last change on this file since 14730 was 14730, checked in by francesca, 3 years ago

lbc_lnk reordering in stp_MLF and other adjustments - ticket #2607

File size: 17.1 KB
Line 
1MODULE stprk3
2   !!======================================================================
3   !!                       ***  MODULE stprk3  ***
4   !! Time-stepping   : manager of the shallow water equation time stepping
5   !!                   3rd order Runge-Kutta scheme
6   !!======================================================================
7   !! History :  NEMO !  2020-03  (A. Nasser, G. Madec)  Original code from  4.0.2
8   !!             -   !  2020-10  (S. Techene, G. Madec)  cleanning
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   stp_RK3       : RK3 Shallow Water Eq. time-stepping
13   !!----------------------------------------------------------------------
14   USE stp_oce        ! modules used in nemo_init and stp_RK3
15   !
16   USE domqco         ! quasi-eulerian coordinate
17   USE phycst         ! physical constants
18   USE usrdef_nam     ! user defined namelist parameters
19
20   IMPLICIT NONE
21   PRIVATE
22
23   PUBLIC   stp_RK3   ! called by nemogcm.F90
24
25   !                                          !**  time level indices  **!
26   INTEGER, PUBLIC ::   Nbb, Nnn, Naa, Nrhs   !: used by nemo_init
27     
28   !! * Substitutions
29#  include "do_loop_substitute.h90"
30#  include "domzgr_substitute.h90"
31   !!----------------------------------------------------------------------
32   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
33   !! $Id: step.F90 12614 2020-03-26 14:59:52Z gm $
34   !! Software governed by the CeCILL license (see ./LICENSE)
35   !!----------------------------------------------------------------------
36CONTAINS
37
38   SUBROUTINE stp_RK3( kstp )
39      !!----------------------------------------------------------------------
40      !!                     ***  ROUTINE stp_RK3  ***
41      !!
42      !! ** Purpose : - RK3 Time stepping scheme for shallow water Eq.
43      !!
44      !! ** Method  : 3rd order time stepping scheme which has 3 stages
45      !!       * Update calendar and forcings
46      !!       stage 1   : n ==> n+1/3 using variables at n
47      !!                 - Compute the rhs of momentum
48      !!                 - Time step ssh at Naa (n+1/3)
49      !!                 - Time step u,v at Naa (n+1/3)
50      !!                 - Swap time indices
51      !!       stage 2   : n ==> n+1/2 using variables at n and n+1/3
52      !!                 - Compute the rhs of momentum
53      !!                 - Time step ssh at Naa (n+1/2)
54      !!                 - Time step u,v at Naa (n+1/2)
55      !!                 - Swap time indices
56      !!       stage 3   : n ==> n+1 using variables at n and n+1/2
57      !!                 - Compute the rhs of momentum
58      !!                 - Time step ssh at Naa (n+1)
59      !!                 - Time step u,v at Naa (n+1)
60      !!                 - Swap time indices
61      !!       * Outputs and diagnostics
62      !!
63      !!       NB: in stages 1 and 2 lateral mixing and forcing are not taken
64      !!           into account in the momentum RHS execpt if key_RK3all is used
65      !!----------------------------------------------------------------------
66      INTEGER, INTENT(in   ) ::   kstp   ! ocean time-step index
67      !
68      INTEGER ::   ji, jj, jk   ! dummy loop indice
69      INTEGER ::   indic        ! error indicator if < 0
70      REAL(wp)::   z1_2rho0,  z5_6,  z3_4  ! local scalars
71      REAL(wp)::   zue3a, zue3b, zua, zrhs_u    ! local scalars
72      REAL(wp)::   zve3a, zve3b, zva, zrhs_v    !   -      -
73      !! ---------------------------------------------------------------------
74      !
75      IF( ln_timing )   CALL timing_start('stp_RK3')
76      !
77      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
78      ! model timestep
79      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
80      !
81      IF ( kstp == nit000 )   ww(:,:,:) = 0._wp   ! initialize vertical velocity one for all to zero
82
83      !
84      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
85      ! update I/O and calendar
86      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
87                             indic = 0                ! reset to no error condition
88                             
89      IF( kstp == nit000 ) THEN                       ! initialize IOM context
90                             CALL iom_init( cxios_context, ld_closedef=.FALSE. )   ! for model grid (including possible AGRIF zoom)
91                             CALL iom_init_closedef
92      ENDIF
93      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
94                             CALL iom_setkt( kstp - nit000 + 1,      cxios_context          )   ! tell IOM we are at time step kstp
95
96      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
97      ! Update external forcing   (SWE: surface boundary condition only)
98      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
99
100                             CALL sbc     ( kstp, Nbb, Nnn )                   ! Sea Boundary Condition
101
102      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
103      ! Ocean physics update   (SWE: eddy viscosity only)
104      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
105
106      IF( l_ldfdyn_time )   CALL ldf_dyn( kstp, Nbb )                          ! eddy viscosity coeff.
107
108      !======================================================================
109      !======================================================================
110      !                     =====       RK3       =====
111      !======================================================================
112      !======================================================================
113
114     
115      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
116      !  RK3 1st stage Ocean dynamics : u, v, ssh
117      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
118      rDt   = rn_Dt / 3._wp 
119      r1_Dt = 1._wp / rDt
120      !
121      !                                 !==  RHS of the momentum Eq.  ==!
122      !
123      uu(:,:,:,Nrhs) = 0._wp                        ! set dynamics trends to zero
124      vv(:,:,:,Nrhs) = 0._wp
125
126      CALL dyn_adv( kstp, Nbb, Nbb, uu, vv, Nrhs )  ! advection (VF or FF) ==> RHS
127      CALL dyn_vor( kstp,      Nbb, uu, vv, Nrhs )  ! vorticity            ==> RHS
128#if defined key_RK3all 
129      CALL dyn_ldf( kstp, Nbb, Nbb, uu, vv, Nrhs )  ! lateral mixing
130#endif
131      z5_6 = 5._wp/6._wp
132      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
133         !                                          ! horizontal pressure gradient
134         zrhs_u =        - grav    * ( ssh(ji+1,jj,Nbb) - ssh(ji,jj,Nbb) ) * r1_e1u(ji,jj)
135         zrhs_v =        - grav    * ( ssh(ji,jj+1,Nbb) - ssh(ji,jj,Nbb) ) * r1_e2v(ji,jj)
136#if defined key_RK3all
137         !                                          ! wind stress and layer friction
138         zrhs_u = zrhs_u + r1_rho0 * ( z5_6*utau_b(ji,jj) + (1._wp - z5_6)*utau(ji,jj) ) / e3u(ji,jj,jk,Nbb)   &
139            &            - rn_rfr  * uu(ji,jj,jk,Nbb)
140         zrhs_v = zrhs_v + r1_rho0 * ( z5_6*vtau_b(ji,jj) + (1._wp - z5_6)*vtau(ji,jj) ) / e3v(ji,jj,jk,Nbb)   &
141            &            - rn_rfr  * vv(ji,jj,jk,Nbb)
142#endif
143         !                                          ! ==> RHS
144         uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + zrhs_u
145         vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + zrhs_v
146      END_3D
147      !
148      !                                 !==  Time stepping of ssh Eq.  ==!   (and update r3_Naa)
149      !
150      CALL ssh_nxt( kstp, Nbb, Nbb, ssh, Naa )      ! after ssh
151      !                                             ! after ssh/h_0 ratio
152      CALL dom_qco_r3c( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) )
153      !
154      !                                 !==  Time stepping of momentum Eq.  ==!
155      !
156      IF( ln_dynadv_vec ) THEN                      ! vector invariant form : applied on velocity
157         DO_3D( 0, 0, 0, 0, 1,jpkm1)
158            uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
159            vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
160         END_3D         
161      ELSE
162         DO_3D( 0, 0, 0, 0, 1,jpkm1)                 ! flux form : applied on thickness weighted velocity
163            zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb)
164            zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb)
165            zue3a = zue3b + rDt * e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
166            zve3a = zve3b + rDt * e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
167            !
168            uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa)   
169            vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa)
170         END_3D
171      ENDIF
172      !
173      CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. )
174      ! [ comm_cleanup ] ! lbc_lnk from DYN - needed for ssh_nxt
175      IF (nn_hls.eq.2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.)
176      !
177      !                                 !==  Swap time levels  ==!
178      Nrhs= Nnn
179      Nnn = Naa
180      Naa = Nrhs
181
182      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
183      !  RK3 2nd stage Ocean dynamics : hdiv, ssh, e3, u, v, w
184      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
185      rDt   = rn_Dt / 2._wp 
186      r1_Dt = 1._wp / rDt
187      !
188      !                                 !==  RHS of the momentum Eq.  ==!
189      !
190      uu(:,:,:,Nrhs) = 0._wp                        ! set dynamics trends to zero
191      vv(:,:,:,Nrhs) = 0._wp
192      CALL dyn_adv( kstp, Nbb, Nnn, uu, vv, Nrhs )  ! advection (VF or FF) ==> RHS
193      CALL dyn_vor( kstp,      Nnn, uu, vv, Nrhs )  ! vorticity            ==> RHS
194#if defined key_RK3all 
195      CALL dyn_ldf( kstp, Nbb, Nbb, uu, vv, Nrhs )  ! lateral mixing
196#endif
197      !
198      z3_4 = 3._wp/4._wp
199      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
200         !                                          ! horizontal pressure gradient
201         zrhs_u =        - grav    * ( ssh(ji+1,jj,Nnn) - ssh(ji,jj,Nnn) ) * r1_e1u(ji,jj)
202         zrhs_v =        - grav    * ( ssh(ji,jj+1,Nnn) - ssh(ji,jj,Nnn) ) * r1_e2v(ji,jj)
203#if defined key_RK3all
204         !                                          ! wind stress and layer friction
205         zrhs_u = zrhs_u + r1_rho0 * ( z3_4*utau_b(ji,jj) + (1._wp - z3_4)*utau(ji,jj) ) / e3u(ji,jj,jk,Nnn)   &
206            &            - rn_rfr  * uu(ji,jj,jk,Nbb)
207         zrhs_v = zrhs_v + r1_rho0 * ( z3_4*vtau_b(ji,jj) + (1._wp - z3_4)*vtau(ji,jj) ) / e3v(ji,jj,jk,Nnn)   &
208            &            - rn_rfr  * vv(ji,jj,jk,Nbb)
209#endif
210         !                                          ! ==> RHS
211         uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + zrhs_u
212         vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + zrhs_v
213      END_3D
214      !
215      !                                 !==  Time stepping of ssh Eq.  ==!   (and update r3_Naa)
216      !
217      CALL ssh_nxt( kstp, Nbb, Nnn, ssh, Naa )      ! after ssh
218      !                                             ! after ssh/h_0 ratio
219      CALL dom_qco_r3c( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) )
220      !
221      !                                 !==  Time stepping of momentum Eq.  ==!
222      !
223      IF( ln_dynadv_vec ) THEN                      ! vector invariant form : applied on velocity
224         DO_3D( 0, 0, 0, 0, 1,jpkm1)
225            uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
226            vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
227         END_3D         
228      ELSE
229         DO_3D( 0, 0, 0, 0, 1,jpkm1)                 ! flux form : applied on thickness weighted velocity
230            zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb)
231            zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb)
232            zue3a = zue3b + rDt * e3u(ji,jj,jk,Nnn) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
233            zve3a = zve3b + rDt * e3v(ji,jj,jk,Nnn) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
234            !
235            uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa)   
236            vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa)
237         END_3D
238      ENDIF
239      !
240      CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. )
241      ! [ comm_cleanup ] ! lbc_lnk from DYN - needed for ssh_nxt
242      IF (nn_hls.eq.2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.)
243      !
244      !                                 !==  Swap time levels  ==!
245      Nrhs= Nnn
246      Nnn = Naa
247      Naa = Nrhs
248       
249      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
250      !  RK3 3rd stage Ocean dynamics : hdiv, ssh, e3, u, v, w
251      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
252      rDt   = rn_Dt
253      r1_Dt = 1._wp / rDt
254      !
255      !                                 !==  RHS of the momentum Eq.  ==!
256      !
257      uu(:,:,:,Nrhs) = 0._wp                        ! set dynamics trends to zero
258      vv(:,:,:,Nrhs) = 0._wp
259      !
260      CALL dyn_adv( kstp, Nbb, Nnn, uu, vv, Nrhs )  ! advection (VF or FF) ==> RHS
261      CALL dyn_vor( kstp,      Nnn, uu, vv, Nrhs )  ! vorticity            ==> RHS
262      CALL dyn_ldf( kstp, Nbb, Nnn, uu, vv, Nrhs )  ! lateral mixing
263
264      z1_2rho0 = 0.5_wp * r1_rho0
265      DO_3D( 0, 0, 0, 0, 1,jpkm1 )
266         !                                          ! horizontal pressure gradient
267         zrhs_u =        - grav * ( ssh(ji+1,jj,Nnn) - ssh(ji,jj,Nnn) ) * r1_e1u(ji,jj)
268         zrhs_v =        - grav * ( ssh(ji,jj+1,Nnn) - ssh(ji,jj,Nnn) ) * r1_e2v(ji,jj)
269         !                                          ! wind stress and layer friction
270         zrhs_u = zrhs_u + z1_2rho0 * ( utau_b(ji,jj) + utau(ji,jj) ) / e3u(ji,jj,jk,Nnn)   &
271            &            - rn_rfr   * uu(ji,jj,jk,Nbb)
272         zrhs_v = zrhs_v + z1_2rho0 * ( vtau_b(ji,jj) + vtau(ji,jj) ) / e3v(ji,jj,jk,Nnn)   &
273            &            - rn_rfr   * vv(ji,jj,jk,Nbb)
274         !                                          ! ==> RHS
275         uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + zrhs_u
276         vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + zrhs_v
277      END_3D
278      !
279      !                                 !==  Time stepping of ssh Eq.  ==!   (and update r3_Naa)
280      !
281      CALL ssh_nxt( kstp, Nbb, Nnn, ssh, Naa )      ! after ssh
282      !                                             ! after ssh/h_0 ratio
283      CALL dom_qco_r3c( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) )
284      !
285      !                                 !==  Time stepping of momentum Eq.  ==!
286      !
287      IF( ln_dynadv_vec ) THEN                      ! vector invariant form : applied on velocity
288         DO_3D( 0, 0, 0, 0, 1,jpkm1)
289            uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
290            vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
291         END_3D
292         !
293      ELSE                                          ! flux form : applied on thickness weighted velocity
294         DO_3D( 0, 0, 0, 0, 1,jpkm1)
295            zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb)
296            zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb)
297            zue3a = zue3b + rDt * e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
298            zve3a = zve3b + rDt * e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
299            !
300            uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa)   
301            vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa)
302         END_3D
303      ENDIF
304      !
305      CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. )
306      ! [ comm_cleanup ] ! lbc_lnk from DYN - needed for ssh_nxt
307      IF (nn_hls.eq.2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.)
308      !
309      !                                 !==  Swap time levels  ==!
310      !
311      Nrhs = Nbb
312      Nbb = Naa
313      Naa = Nrhs
314
315      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
316      ! diagnostics and outputs at Nbb (i.e. the just computed time step)
317      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
318     
319      IF( ln_diacfl  )   CALL dia_cfl      ( kstp,      Nbb )      ! Courant number diagnostics
320                         CALL dia_wri      ( kstp,      Nbb )      ! ocean model: outputs
321      !
322      IF( lrst_oce   )   CALL rst_write    ( kstp, Nbb, Nbb )   ! write output ocean restart file
323
324      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
325      ! Control
326      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
327                         CALL stp_ctl      ( kstp     , Nbb )
328
329      IF( kstp == nit000 ) THEN                          ! 1st time step only
330                                        CALL iom_close( numror )   ! close input  ocean restart file
331         IF(lwm)                        CALL FLUSH    ( numond )   ! flush output namelist oce
332         IF(lwm .AND. numoni /= -1 )    CALL FLUSH    ( numoni )   ! flush output namelist ice (if exist)
333      ENDIF
334
335      !
336#if defined key_xios
337      IF( kstp == nitend .OR. indic < 0 ) THEN
338         CALL iom_context_finalize( cxios_context )
339      ENDIF
340#endif
341      !
342      IF( ln_timing )   CALL timing_stop('stp_RK3')
343      !
344   END SUBROUTINE stp_RK3
345
346   !!======================================================================
347END MODULE stprk3
Note: See TracBrowser for help on using the repository browser.