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_MERGE_2013/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 4317

Last change on this file since 4317 was 4313, checked in by cetlod, 10 years ago

dev_MERGE_2013 : fix on argument call of eos

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