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.
step.F90 in NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/SWE – NEMO

source: NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/SWE/step.F90 @ 13553

Last change on this file since 13553 was 13553, checked in by hadcv, 3 years ago

Merge in trunk up to [13550]

File size: 16.8 KB
Line 
1MODULE step
2   !!======================================================================
3   !!                       ***  MODULE step  ***
4   !! Time-stepping   : manager of the shallow water equation time stepping
5   !!======================================================================
6   !! History :  NEMO !  2020-03  (A. Nasser, G. Madec)  Original code from  4.0.2
7   !!----------------------------------------------------------------------
8#if defined key_qco
9   !!----------------------------------------------------------------------
10   !!   'key_qco'      EMPTY MODULE      Quasi-Eulerian vertical coordonate
11   !!----------------------------------------------------------------------
12#else
13   !!----------------------------------------------------------------------
14   !!   stp             : Shallow Water time-stepping
15   !!----------------------------------------------------------------------
16   USE step_oce         ! time stepping definition modules
17   USE phycst           ! physical constants
18   USE usrdef_nam
19   !
20   USE iom              ! xIOs server
21
22   IMPLICIT NONE
23   PRIVATE
24
25   PUBLIC   stp   ! called by nemogcm.F90
26   
27   !!----------------------------------------------------------------------
28   !! time level indices
29   !!----------------------------------------------------------------------
30   INTEGER, PUBLIC ::   Nbb, Nnn, Naa, Nrhs   !! used by nemo_init
31     
32   !! * Substitutions
33#  include "do_loop_substitute.h90"
34   !!----------------------------------------------------------------------
35   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
36   !! $Id: step.F90 12614 2020-03-26 14:59:52Z gm $
37   !! Software governed by the CeCILL license (see ./LICENSE)
38   !!----------------------------------------------------------------------
39CONTAINS
40
41#if defined key_agrif
42   RECURSIVE SUBROUTINE stp( )
43      INTEGER             ::   kstp   ! ocean time-step index
44#else
45   SUBROUTINE stp( kstp )
46      INTEGER, INTENT(in) ::   kstp   ! ocean time-step index
47#endif
48      !!----------------------------------------------------------------------
49      !!                     ***  ROUTINE stp  ***
50      !!
51      !! ** Purpose : - Time stepping of shallow water (SHW) (momentum and ssh eqs.)
52      !!
53      !! ** Method  : -1- Update forcings
54      !!              -2- Update the ssh at Naa
55      !!              -3- Compute the momentum trends (Nrhs)
56      !!              -4- Update the horizontal velocity
57      !!              -5- Apply Asselin time filter to uu,vv,ssh
58      !!              -6- Outputs and diagnostics
59      !!----------------------------------------------------------------------
60      INTEGER ::   ji, jj, jk   ! dummy loop indice
61      INTEGER ::   indic        ! error indicator if < 0
62!!gm kcall can be removed, I guess
63      INTEGER ::   kcall        ! optional integer argument (dom_vvl_sf_nxt)
64      REAL(wp)::   z1_2rho0     ! local scalars
65     
66      REAL(wp) ::   zue3a, zue3n, zue3b    ! local scalars
67      REAL(wp) ::   zve3a, zve3n, zve3b    !   -      -
68      REAL(wp) ::   ze3t_tf, ze3u_tf, ze3v_tf, zua, zva
69      !! ---------------------------------------------------------------------
70#if defined key_agrif
71      kstp = nit000 + Agrif_Nb_Step()
72      Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs   ! agrif_oce module copies of time level indices
73      IF( lk_agrif_debug ) THEN
74         IF( Agrif_Root() .and. lwp)   WRITE(*,*) '---'
75         IF(lwp)   WRITE(*,*) 'Grid Number', Agrif_Fixed(),' time step ', kstp, 'int tstep', Agrif_NbStepint()
76      ENDIF
77      IF( kstp == nit000 + 1 )   lk_agrif_fstep = .FALSE.
78# if defined key_iomput
79      IF( Agrif_Nbstepint() == 0 )   CALL iom_swap( cxios_context )
80# endif
81#endif
82      !
83      IF( ln_timing )   CALL timing_start('stp')
84      !
85      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
86      ! model timestep
87      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
88      !
89      IF( l_1st_euler ) THEN     ! start or restart with Euler 1st time-step
90         rDt   =  rn_Dt   
91         r1_Dt = 1._wp / rDt
92      ENDIF
93
94      IF ( kstp == nit000 )   ww(:,:,:) = 0._wp   ! initialize vertical velocity one for all to zero
95
96      !
97      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
98      ! update I/O and calendar
99      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
100                             indic = 0                ! reset to no error condition
101                             
102      IF( kstp == nit000 ) THEN                       ! initialize IOM context (must be done after nemo_init for AGRIF+XIOS+OASIS)
103                             CALL iom_init( cxios_context, ld_closedef=.FALSE. )   ! for model grid (including passible AGRIF zoom)
104         IF( lk_diamlr   )   CALL dia_mlr_iom_init    ! with additional setup for multiple-linear-regression analysis
105                             CALL iom_init_closedef
106         IF( ln_crs      )   CALL iom_init( TRIM(cxios_context)//"_crs" )  ! for coarse grid
107      ENDIF
108      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
109                             CALL iom_setkt( kstp - nit000 + 1,      cxios_context          )   ! tell IOM we are at time step kstp
110      IF( ln_crs         )   CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" )   ! tell IOM we are at time step kstp
111
112      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
113      ! Update external forcing (tides, open boundaries, ice shelf interaction and surface boundary condition (including sea-ice)
114      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
115      IF( ln_tide    )   CALL tide_update( kstp )                     ! update tide potential
116      IF( ln_apr_dyn )   CALL sbc_apr ( kstp )                        ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib)
117      IF( ln_bdy     )   CALL bdy_dta ( kstp, Nnn )                   ! update dynamic & tracer data at open boundaries
118                         CALL sbc     ( kstp, Nbb, Nnn )              ! Sea Boundary Condition
119
120      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
121      ! Ocean physics update
122      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
123      !  LATERAL  PHYSICS
124      !                                                                        ! eddy diffusivity coeff.
125      IF( l_ldfdyn_time )   CALL ldf_dyn( kstp, Nbb )                          ! eddy viscosity coeff.
126
127      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
128      !  Ocean dynamics : hdiv, ssh, e3, u, v, w
129      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
130
131                            CALL ssh_nxt       ( kstp, Nbb, Nnn, ssh, Naa )    ! after ssh (includes call to div_hor)
132                         uu(:,:,:,Nrhs) = 0._wp            ! set dynamics trends to zero
133                         vv(:,:,:,Nrhs) = 0._wp
134
135      IF( .NOT.ln_linssh )  CALL dom_vvl_sf_nxt( kstp, Nbb, Nnn,      Naa )    ! after vertical scale factors
136
137      IF( ln_bdy     )      CALL bdy_dyn3d_dmp ( kstp, Nbb,      uu, vv, Nrhs )  ! bdy damping trends
138
139#if defined key_agrif
140      IF(.NOT. Agrif_Root())  & 
141               &            CALL Agrif_Sponge_dyn        ! momentum sponge
142#endif
143                            CALL dyn_adv( kstp, Nbb, Nnn      , uu, vv, Nrhs )  ! advection (VF or FF)   ==> RHS
144 
145                            CALL dyn_vor( kstp,      Nnn      , uu, vv, Nrhs )  ! vorticity              ==> RHS
146 
147                            CALL dyn_ldf( kstp, Nbb, Nnn      , uu, vv, Nrhs )  ! lateral mixing
148
149!!an - calcul du gradient de pression horizontal (explicit)
150      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
151         uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) - grav * ( ssh(ji+1,jj,Nnn) - ssh(ji,jj,Nnn) ) * r1_e1u(ji,jj)
152         vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) - grav * ( ssh(ji,jj+1,Nnn) - ssh(ji,jj,Nnn) ) * r1_e2v(ji,jj)
153      END_3D
154      !
155      ! add wind stress forcing and layer linear friction to the RHS
156      z1_2rho0 = 0.5_wp * r1_rho0
157      DO_3D( 0, 0, 0, 0,1,jpkm1)
158         uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + z1_2rho0 * ( utau_b(ji,jj) + utau(ji,jj) ) / e3u(ji,jj,jk,Nnn)   &
159            &                                  - rn_rfr * uu(ji,jj,jk,Nbb)
160         vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + z1_2rho0 * ( vtau_b(ji,jj) + vtau(ji,jj) ) / e3v(ji,jj,jk,Nnn)   &
161            &                                  - rn_rfr * vv(ji,jj,jk,Nbb)
162      END_3D   
163!!an         
164 
165      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
166      ! Leap-Frog time splitting + Robert-Asselin time filter on u,v,e3
167      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
168         
169!! what about  IF( .NOT.ln_linssh )  ?
170!!an futur module dyn_nxt (a la place de dyn_atf)
171     
172      IF( ln_dynadv_vec ) THEN      ! vector invariant form : applied on velocity
173         IF( l_1st_euler ) THEN           ! Euler time stepping (no Asselin filter)
174            DO_3D( 0, 0, 0, 0,1,jpkm1)
175               uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
176               vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
177             END_3D         
178         ELSE                             ! Leap Frog time stepping + Asselin filter         
179            DO_3D( 1, 1, 1, 1,1,jpkm1)
180               zua = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
181               zva = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
182               !                                                                  ! Asselin time filter on u,v (Nnn)
183               uu(ji,jj,jk,Nnn) = uu(ji,jj,jk,Nnn) + rn_atfp * (uu(ji,jj,jk,Nbb) - 2._wp * uu(ji,jj,jk,Nnn) + zua)
184               vv(ji,jj,jk,Nnn) = vv(ji,jj,jk,Nnn) + rn_atfp * (vv(ji,jj,jk,Nbb) - 2._wp * vv(ji,jj,jk,Nnn) + zva)
185               !
186               ze3u_tf = e3u(ji,jj,jk,Nnn) + rn_atfp * ( e3u(ji,jj,jk,Nbb) - 2._wp * e3u(ji,jj,jk,Nnn)  + e3u(ji,jj,jk,Naa) )
187               ze3v_tf = e3v(ji,jj,jk,Nnn) + rn_atfp * ( e3v(ji,jj,jk,Nbb) - 2._wp * e3v(ji,jj,jk,Nnn)  + e3v(ji,jj,jk,Naa) )
188               ze3t_tf = e3t(ji,jj,jk,Nnn) + rn_atfp * ( e3t(ji,jj,jk,Nbb) - 2._wp * e3t(ji,jj,jk,Nnn)  + e3t(ji,jj,jk,Naa) )
189               !
190               e3u(ji,jj,jk,Nnn) = ze3u_tf
191               e3v(ji,jj,jk,Nnn) = ze3v_tf
192               e3t(ji,jj,jk,Nnn) = ze3t_tf
193               !             
194               uu(ji,jj,jk,Naa) = zua
195               vv(ji,jj,jk,Naa) = zva
196            END_3D
197         ENDIF
198         !
199      ELSE                          ! flux form : applied on thickness weighted velocity
200         IF( l_1st_euler ) THEN           ! Euler time stepping (no Asselin filter)
201            DO_3D( 0, 0, 0, 0,1,jpkm1)
202               zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb)
203               zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb)
204               !                                                ! LF time stepping
205               zue3a = zue3b + rDt * e3u(ji,jj,jk,Nrhs) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
206               zve3a = zve3b + rDt * e3v(ji,jj,jk,Nrhs) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
207               !
208               uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa)   
209               vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa)
210            END_3D
211         ELSE                             ! Leap Frog time stepping + Asselin filter
212            DO_3D( 1, 1, 1, 1,1,jpkm1)
213               zue3n = e3u(ji,jj,jk,Nnn) * uu(ji,jj,jk,Nnn)
214               zve3n = e3v(ji,jj,jk,Nnn) * vv(ji,jj,jk,Nnn)
215               zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb)
216               zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb)
217               !                                                ! LF time stepping
218               zue3a = zue3b + rDt * e3u(ji,jj,jk,Nrhs) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk)
219               zve3a = zve3b + rDt * e3v(ji,jj,jk,Nrhs) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk)
220               !                                                ! Asselin time filter on e3u/v/t
221               ze3u_tf = e3u(ji,jj,jk,Nnn) + rn_atfp * ( e3u(ji,jj,jk,Nbb) - 2._wp * e3u(ji,jj,jk,Nnn)  + e3u(ji,jj,jk,Naa) ) 
222               ze3v_tf = e3v(ji,jj,jk,Nnn) + rn_atfp * ( e3v(ji,jj,jk,Nbb) - 2._wp * e3v(ji,jj,jk,Nnn)  + e3v(ji,jj,jk,Naa) )
223               ze3t_tf = e3t(ji,jj,jk,Nnn) + rn_atfp * ( e3t(ji,jj,jk,Nbb) - 2._wp * e3t(ji,jj,jk,Nnn)  + e3t(ji,jj,jk,Naa) ) 
224               !                                                ! Asselin time filter on u,v (Nnn)
225               uu(ji,jj,jk,Nnn) = ( zue3n + rn_atfp * ( zue3b - 2._wp * zue3n  + zue3a ) ) / ze3u_tf
226               vv(ji,jj,jk,Nnn) = ( zve3n + rn_atfp * ( zve3b - 2._wp * zve3n  + zve3a ) ) / ze3v_tf
227               !
228               e3u(ji,jj,jk,Nnn) = ze3u_tf
229               e3v(ji,jj,jk,Nnn) = ze3v_tf
230               e3t(ji,jj,jk,Nnn) = ze3t_tf
231               !
232               uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa)   
233               vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa)
234            END_3D       
235         ENDIF
236      ENDIF
237
238
239      CALL lbc_lnk_multi( 'stp', uu(:,:,:,Nnn), 'U', -1., vv(:,:,:,Nnn), 'V', -1.,   &   !* local domain boundaries
240         &                       uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1.    )     
241
242!!an         
243
244      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
245      ! Set boundary conditions, time filter and swap time levels
246      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
247
248!!an TO BE ADDED : dyn_nxt
249!!                         CALL dyn_atf       ( kstp, Nbb, Nnn, Naa, uu, vv, e3t, e3u, e3v  )  ! time filtering of "now" velocities and scale factors
250!!an TO BE ADDED : a simplifier
251!!                         CALL ssh_atf       ( kstp, Nbb, Nnn, Naa, ssh )                     ! time filtering of "now" sea surface height
252      IF ( .NOT.( l_1st_euler ) ) THEN   ! Only do time filtering for leapfrog timesteps
253         !                                                  ! filtering "now" field
254         ssh(:,:,Nnn) = ssh(:,:,Nnn) + rn_atfp * ( ssh(:,:,Nbb) - 2 * ssh(:,:,Nnn) + ssh(:,:,Naa) )
255      ENDIF
256!!an
257
258
259      ! Swap time levels
260      Nrhs = Nbb
261      Nbb = Nnn
262      Nnn = Naa
263      Naa = Nrhs
264      !
265                         CALL dom_vvl_sf_update( kstp, Nbb, Nnn, Naa )  ! recompute vertical scale factors
266      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
267      ! diagnostics and outputs
268      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
269      IF( ln_floats  )   CALL flo_stp   ( kstp, Nbb, Nnn )      ! drifting Floats
270      IF( ln_diacfl  )   CALL dia_cfl   ( kstp,      Nnn )      ! Courant number diagnostics
271   
272                         CALL dia_wri   ( kstp,      Nnn )      ! ocean model: outputs
273      !
274      IF( lrst_oce   )   CALL rst_write    ( kstp, Nbb, Nnn )   ! write output ocean restart file
275         
276
277#if defined key_agrif
278      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
279      ! AGRIF
280      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<     
281                         Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs      ! agrif_oce module copies of time level indices
282                         CALL Agrif_Integrate_ChildGrids( stp )       ! allows to finish all the Child Grids before updating
283
284                         IF( Agrif_NbStepint() == 0 ) THEN
285                            CALL Agrif_update_all( )                  ! Update all components
286                         ENDIF
287#endif
288      IF( ln_diaobs  )   CALL dia_obs      ( kstp, Nnn )      ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
289
290      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
291      ! Control
292      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
293                         CALL stp_ctl      ( kstp, Nnn )
294               
295                         
296      IF( kstp == nit000 ) THEN                          ! 1st time step only
297                                        CALL iom_close( numror )   ! close input  ocean restart file
298         IF(lwm)                        CALL FLUSH    ( numond )   ! flush output namelist oce
299         IF(lwm .AND. numoni /= -1 )    CALL FLUSH    ( numoni )   ! flush output namelist ice (if exist)
300      ENDIF
301
302      !
303#if defined key_iomput
304      IF( kstp == nitend .OR. indic < 0 ) THEN
305                      CALL iom_context_finalize(      cxios_context          ) ! needed for XIOS+AGRIF
306                      IF(lrxios) CALL iom_context_finalize(      crxios_context          )
307      ENDIF
308#endif
309      !
310      IF( l_1st_euler ) THEN         ! recover Leap-frog timestep
311         rDt = 2._wp * rn_Dt   
312         r1_Dt = 1._wp / rDt
313         l_1st_euler = .FALSE.     
314      ENDIF
315      !
316      IF( ln_timing )   CALL timing_stop('stp')
317      !
318   END SUBROUTINE stp
319#endif
320   !
321   !!======================================================================
322END MODULE step
Note: See TracBrowser for help on using the repository browser.