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 branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 4267

Last change on this file since 4267 was 4258, checked in by acc, 10 years ago

Branch 2013/dev_r3858_NOC_ZTC, #863. Merge in trunk changes from 3858 to 4119; resolve conflicts and check in prior to compiling and testing

  • Property svn:keywords set to Id
File size: 19.5 KB
Line 
1MODULE step
2   !!======================================================================
3   !!                       ***  MODULE step  ***
4   !! Time-stepping    : manager of the ocean, tracer and ice time stepping
5   !!======================================================================
6   !! History :  OPA  !  1991-03  (G. Madec)  Original code
7   !!             -   !  1991-11  (G. Madec)
8   !!             -   !  1992-06  (M. Imbard)  add a first output record
9   !!             -   !  1996-04  (G. Madec)  introduction of dynspg
10   !!             -   !  1996-04  (M.A. Foujols)  introduction of passive tracer
11   !!            8.0  !  1997-06  (G. Madec)  new architecture of call
12   !!            8.2  !  1997-06  (G. Madec, M. Imbard, G. Roullet)  free surface
13   !!             -   !  1999-02  (G. Madec, N. Grima)  hpg implicit
14   !!             -   !  2000-07  (J-M Molines, M. Imbard)  Open Bondary Conditions
15   !!   NEMO     1.0  !  2002-06  (G. Madec)  free form, suppress macro-tasking
16   !!             -   !  2004-08  (C. Talandier) New trends organization
17   !!             -   !  2005-01  (C. Ethe) Add the KPP closure scheme
18   !!             -   !  2005-11  (G. Madec)  Reorganisation of tra and dyn calls
19   !!             -   !  2006-01  (L. Debreu, C. Mazauric)  Agrif implementation
20   !!             -   !  2006-07  (S. Masson)  restart using iom
21   !!            3.2  !  2009-02  (G. Madec, R. Benshila)  reintroduicing z*-coordinate
22   !!             -   !  2009-06  (S. Masson, G. Madec)  TKE restart compatible with key_cpl
23   !!            3.3  !  2010-05  (K. Mogensen, A. Weaver, M. Martin, D. Lea) Assimilation interface
24   !!             -   !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA
25   !!            3.4  !  2011-04  (G. Madec, C. Ethe) Merge of dtatem and dtasal
26   !!----------------------------------------------------------------------
27
28   !!----------------------------------------------------------------------
29   !!   stp             : OPA system time-stepping
30   !!----------------------------------------------------------------------
31   USE step_oce         ! time stepping definition modules
32
33   IMPLICIT NONE
34   PRIVATE
35
36   PUBLIC   stp   ! called by opa.F90
37
38   !! * Substitutions
39#  include "domzgr_substitute.h90"
40#  include "zdfddm_substitute.h90"
41   !!----------------------------------------------------------------------
42   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
43   !! $Id$
44   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
45   !!----------------------------------------------------------------------
46CONTAINS
47
48#if defined key_agrif
49   SUBROUTINE stp( )
50      INTEGER             ::   kstp   ! ocean time-step index
51#else
52   SUBROUTINE stp( kstp )
53      INTEGER, INTENT(in) ::   kstp   ! ocean time-step index
54#endif
55      !!----------------------------------------------------------------------
56      !!                     ***  ROUTINE stp  ***
57      !!
58      !! ** Purpose : - Time stepping of OPA (momentum and active tracer eqs.)
59      !!              - Time stepping of LIM (dynamic and thermodynamic eqs.)
60      !!              - Tme stepping  of TRC (passive tracer eqs.)
61      !!
62      !! ** Method  : -1- Update forcings and data
63      !!              -2- Update ocean physics
64      !!              -3- Compute the t and s trends
65      !!              -4- Update t and s
66      !!              -5- Compute the momentum trends
67      !!              -6- Update the horizontal velocity
68      !!              -7- Compute the diagnostics variables (rd,N2, div,cur,w)
69      !!              -8- Outputs and diagnostics
70      !!----------------------------------------------------------------------
71      INTEGER ::   jk       ! dummy loop indice
72      INTEGER ::   indic    ! error indicator if < 0
73      !! ---------------------------------------------------------------------
74
75#if defined key_agrif
76      kstp = nit000 + Agrif_Nb_Step()
77!      IF ( Agrif_Root() .and. lwp) Write(*,*) '---'
78!      IF (lwp) Write(*,*) 'Grid Number',Agrif_Fixed(),' time step ',kstp
79# if defined key_iomput
80      IF( Agrif_Nbstepint() == 0 )   CALL iom_swap
81# endif
82#endif
83                             indic = 0                ! reset to no error condition
84      IF( kstp == nit000 )   CALL iom_init            ! iom_put initialization (must be done after nemo_init for AGRIF+XIOS+OASIS)
85      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
86                             CALL iom_setkt( kstp - nit000 + 1 )   ! say to iom that we are at time step kstp
87
88      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
89      ! Update data, open boundaries, surface boundary condition (including sea-ice)
90      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
91                         CALL sbc     ( kstp )        ! Sea Boundary Condition (including sea-ice)
92      IF( lk_tide    )   CALL sbc_tide( kstp )
93      IF( lk_obc     )   CALL obc_dta ( kstp )        ! update dynamic and tracer data at open boundaries
94      IF( lk_obc     )   CALL obc_rad ( kstp )        ! compute phase velocities at open boundaries
95      IF( lk_bdy     )   CALL bdy_dta ( kstp, time_offset=+1 )   ! update dynamic & tracer data at open boundaries
96
97      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
98      !  Ocean dynamics : hdiv, rot, ssh, e3, wn
99      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
100                         CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_cur)
101      IF( lk_dynspg_ts ) THEN
102                                  CALL wzv           ( kstp )  ! now cross-level velocity
103          ! In case the time splitting case, update almost all momentum trends here:
104          ! Note that the computation of vertical velocity above, hence "after" sea level
105          ! is necessary to compute momentum advection for the rhs of barotropic loop:
106                                  CALL eos    ( tsn, rhd, rhop )                 ! now in situ density for hpg computation
107          IF( ln_zps      )       CALL zps_hde( kstp, jpts, tsn, gtsu, gtsv,  &  ! zps: now hor. derivative
108                &                                          rhd, gru , grv  )     ! of t, s, rd at the last ocean level
109                                  CALL zdf_bfr( kstp )         ! bottom friction (if quadratic)
110
111                                  ua(:,:,:) = 0.e0             ! set dynamics trends to zero
112                                  va(:,:,:) = 0.e0
113          IF(  ln_asmiau .AND. &
114             & ln_dyninc       )  CALL dyn_asm_inc  ( kstp )   ! apply dynamics assimilation increment
115          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! subtract Neptune velocities (simplified)
116          IF( lk_bdy           )  CALL bdy_dyn3d_dmp( kstp )   ! bdy damping trends
117                                  CALL dyn_adv      ( kstp )   ! advection (vector or flux form)
118                                  CALL dyn_vor      ( kstp )   ! vorticity term including Coriolis
119                                  CALL dyn_ldf      ( kstp )   ! lateral mixing
120          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! add Neptune velocities (simplified)
121#if defined key_agrif
122          IF(.NOT. Agrif_Root())  CALL Agrif_Sponge_dyn        ! momentum sponge
123#endif
124                                  CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
125                                  CALL dyn_spg( kstp, indic )  ! surface pressure gradient
126                                  CALL div_cur( kstp )         ! Horizontal divergence & Relative vorticity (2nd call in time-split case)
127
128                                  ua_sv(:,:,:) = ua(:,:,:)     ! save next velocities (not trends !)
129                                  va_sv(:,:,:) = va(:,:,:)
130      ENDIF
131      IF( lk_vvl     )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors
132                         CALL wzv           ( kstp )  ! now cross-level velocity (original)
133
134      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
135      ! Ocean physics update                (ua, va, tsa used as workspace)
136      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
137                         CALL bn2( tsb, rn2b )        ! before Brunt-Vaisala frequency
138                         CALL bn2( tsn, rn2  )        ! now    Brunt-Vaisala frequency
139      !
140      !  VERTICAL PHYSICS
141      IF( .NOT. lk_dynspg_ts ) CALL zdf_bfr( kstp )         ! bottom friction
142
143      !                                               ! Vertical eddy viscosity and diffusivity coefficients
144      IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz
145      IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz
146      IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz
147      IF( lk_zdfkpp  )   CALL zdf_kpp( kstp )            ! KPP closure scheme for Kz
148      IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value)
149         avt (:,:,:) = rn_avt0 * tmask(:,:,:)
150         avmu(:,:,:) = rn_avm0 * umask(:,:,:)
151         avmv(:,:,:) = rn_avm0 * vmask(:,:,:)
152      ENDIF
153      IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths
154         DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:) * tmask(:,:,jk)   ;   END DO
155      ENDIF
156      IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity
157
158      IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing
159
160      IF( lk_zdfddm .AND. .NOT. lk_zdfkpp )   &
161         &               CALL zdf_ddm( kstp )         ! double diffusive mixing
162
163                         CALL zdf_mxl( kstp )         ! mixed layer depth
164
165                                                      ! write TKE or GLS information in the restart file
166      IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' )
167      IF( lrst_oce .AND. lk_zdfgls )   CALL gls_rst( kstp, 'WRITE' )
168      !
169      !  LATERAL  PHYSICS
170      !
171      IF( lk_ldfslp ) THEN                            ! slope of lateral mixing
172                         CALL eos( tsb, rhd, gdept_0(:,:,:) )             ! before in situ density
173         IF( ln_zps )    CALL zps_hde( kstp, jpts, tsb, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
174            &                                      rhd, gru , grv  )      ! of t, s, rd at the last ocean level
175         IF( ln_traldf_grif ) THEN                           ! before slope for Griffies operator
176                         CALL ldf_slp_grif( kstp )
177         ELSE
178                         CALL ldf_slp( kstp, rhd, rn2b )     ! before slope for Madec operator
179         ENDIF
180      ENDIF
181#if defined key_traldf_c2d
182      IF( lk_traldf_eiv )   CALL ldf_eiv( kstp )      ! eddy induced velocity coefficient
183#endif
184#if defined key_traldf_c3d && key_traldf_smag
185                          CALL ldf_tra_smag( kstp )      ! eddy induced velocity coefficient
186#  endif
187#if defined key_dynldf_c3d && key_dynldf_smag
188                          CALL ldf_dyn_smag( kstp )      ! eddy induced velocity coefficient
189#  endif
190
191
192      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
193      ! diagnostics and outputs             (ua, va, tsa used as workspace)
194      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
195      IF( lk_floats  )   CALL flo_stp( kstp )         ! drifting Floats
196      IF( lk_diahth  )   CALL dia_hth( kstp )         ! Thermocline depth (20 degres isotherm depth)
197      IF( lk_diafwb  )   CALL dia_fwb( kstp )         ! Fresh water budget diagnostics
198      IF( ln_diaptr  )   CALL dia_ptr( kstp )         ! Poleward TRansports diagnostics
199      IF( lk_diadct  )   CALL dia_dct( kstp )         ! Transports
200      IF( lk_diaar5  )   CALL dia_ar5( kstp )         ! ar5 diag
201      IF( lk_diaharm )   CALL dia_harm( kstp )        ! Tidal harmonic analysis
202                         CALL dia_wri( kstp )         ! ocean model: outputs
203
204#if defined key_top
205      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
206      ! Passive Tracer Model
207      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
208                         CALL trc_stp( kstp )         ! time-stepping
209#endif
210
211      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
212      ! Active tracers                              (ua, va used as workspace)
213      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
214                             tsa(:,:,:,:) = 0.e0            ! set tracer trends to zero
215
216      IF(  ln_asmiau .AND. &
217         & ln_trainc     )   CALL tra_asm_inc( kstp )       ! apply tracer assimilation increment
218                             CALL tra_sbc    ( kstp )       ! surface boundary condition
219      IF( ln_traqsr      )   CALL tra_qsr    ( kstp )       ! penetrative solar radiation qsr
220      IF( ln_trabbc      )   CALL tra_bbc    ( kstp )       ! bottom heat flux
221      IF( lk_trabbl      )   CALL tra_bbl    ( kstp )       ! advective (and/or diffusive) bottom boundary layer scheme
222      IF( ln_tradmp      )   CALL tra_dmp    ( kstp )       ! internal damping trends
223      IF( lk_bdy         )   CALL bdy_tra_dmp( kstp )       ! bdy damping trends
224                             CALL tra_adv    ( kstp )       ! horizontal & vertical advection
225      IF( lk_zdfkpp      )   CALL tra_kpp    ( kstp )       ! KPP non-local tracer fluxes
226                             CALL tra_ldf    ( kstp )       ! lateral mixing
227#if defined key_agrif
228      IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_tra          ! tracers sponge
229#endif
230                             CALL tra_zdf    ( kstp )       ! vertical mixing and after tracer fields
231
232      IF( ln_dynhpg_imp  ) THEN                             ! semi-implicit hpg (time stepping then eos)
233         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
234                             CALL tra_nxt( kstp )                ! tracer fields at next time step
235                             CALL eos    ( tsa, rhd, rhop, fsdept_n(:,:,:) )  ! Time-filtered in situ density for hpg computation
236         IF( ln_zps      )   CALL zps_hde( kstp, jpts, tsa, gtsu, gtsv,  &    ! zps: time filtered hor. derivative
237            &                                          rhd, gru , grv  )      ! of t, s, rd at the last ocean level
238
239      ELSE                                                  ! centered hpg  (eos then time stepping)
240         IF ( .NOT. lk_dynspg_ts ) THEN                     ! eos already called in time-split case
241                                CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) )  ! now in situ density for hpg computation
242            IF( ln_zps      )   CALL zps_hde( kstp, jpts, tsn, gtsu, gtsv,  &    ! zps: now hor. derivative
243            &                                          rhd, gru , grv  )      ! of t, s, rd at the last ocean level
244         ENDIF
245         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
246                             CALL tra_nxt( kstp )                ! tracer fields at next time step
247      ENDIF
248
249      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
250      ! Dynamics                                    (tsa used as workspace)
251      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
252      IF( lk_dynspg_ts   )  THEN
253! revert to previously computed tendencies:
254! (not using ua, va as temporary arrays during tracers' update could avoid that)
255                               ua(:,:,:) = ua_sv(:,:,:)
256                               va(:,:,:) = va_sv(:,:,:)
257                               CALL dyn_bfr( kstp )         ! bottom friction
258                               CALL dyn_zdf( kstp )         ! vertical diffusion
259      ELSE
260                               ua(:,:,:) = 0.e0             ! set dynamics trends to zero
261                               va(:,:,:) = 0.e0
262
263        IF(  ln_asmiau .AND. &
264           & ln_dyninc      )  CALL dyn_asm_inc( kstp )     ! apply dynamics assimilation increment
265        IF( ln_bkgwri )        CALL asm_bkg_wri( kstp )     ! output background fields
266        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! subtract Neptune velocities (simplified)
267        IF( lk_bdy          )  CALL bdy_dyn3d_dmp(kstp )    ! bdy damping trends
268                               CALL dyn_adv( kstp )         ! advection (vector or flux form)
269                               CALL dyn_vor( kstp )         ! vorticity term including Coriolis
270                               CALL dyn_ldf( kstp )         ! lateral mixing
271        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! add Neptune velocities (simplified)
272#if defined key_agrif
273        IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_dyn        ! momemtum sponge
274#endif
275                               CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
276                               CALL dyn_bfr( kstp )         ! bottom friction
277                               CALL dyn_zdf( kstp )         ! vertical diffusion
278                               CALL dyn_spg( kstp, indic )  ! surface pressure gradient
279      ENDIF
280                               CALL dyn_nxt( kstp )         ! lateral velocity at next time step
281
282                               CALL ssh_swp( kstp )         ! swap of sea surface height
283      IF( lk_vvl           )   CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors
284
285      IF( ln_diahsb        )   CALL dia_hsb( kstp )         ! - ML - global conservation diagnostics
286      IF( lk_diaobs  )         CALL dia_obs( kstp )         ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
287
288      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
289      ! Control and restarts
290      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
291                               CALL stp_ctl( kstp, indic )
292      IF( indic < 0        )   THEN
293                               CALL ctl_stop( 'step: indic < 0' )
294                               CALL dia_wri_state( 'output.abort', kstp )
295      ENDIF
296      IF( kstp == nit000   )   CALL iom_close( numror )     ! close input  ocean restart file
297      IF( lrst_oce         )   CALL rst_write    ( kstp )   ! write output ocean restart file
298      IF( lk_obc           )   CALL obc_rst_write( kstp )   ! write open boundary restart file
299
300      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
301      ! Trends                              (ua, va, tsa used as workspace)
302      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
303      IF( nstop == 0 ) THEN
304         IF( lk_trddyn     )   CALL trd_dwr( kstp )         ! trends: dynamics
305         IF( lk_trdtra     )   CALL trd_twr( kstp )         ! trends: active tracers
306         IF( lk_trdmld     )   CALL trd_mld( kstp )         ! trends: Mixed-layer
307         IF( lk_trdvor     )   CALL trd_vor( kstp )         ! trends: vorticity budget
308      ENDIF
309
310      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
311      ! Coupled mode
312      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
313      IF( lk_cpl           )   CALL sbc_cpl_snd( kstp )     ! coupled mode : field exchanges
314      !
315#if defined key_iomput
316      IF( kstp == nitend .OR. indic < 0 )   CALL xios_context_finalize() ! needed for XIOS+AGRIF
317#endif
318      !
319      IF( nn_timing == 1 .AND.  kstp == nit000  )   CALL timing_reset
320      !
321   END SUBROUTINE stp
322
323   !!======================================================================
324END MODULE step
Note: See TracBrowser for help on using the repository browser.