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/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE – NEMO

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/step.F90 @ 10946

Last change on this file since 10946 was 10946, checked in by acc, 5 years ago

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps : Convert STO, TRD and USR modules and all knock on effects of these conversions. Note change to USR module may have implications for the TEST CASES (not tested yet). Standard SETTE tested only

  • Property svn:keywords set to Id
File size: 22.2 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   !!            3.6  !  2012-07  (J. Simeon, G. Madec. C. Ethe)  Online coarsening of outputs
27   !!            3.6  !  2014-04  (F. Roquet, G. Madec) New equations of state
28   !!            3.6  !  2014-10  (E. Clementi, P. Oddo) Add Qiao vertical mixing in case of waves
29   !!            3.7  !  2014-10  (G. Madec)  LDF simplication
30   !!             -   !  2014-12  (G. Madec) remove KPP scheme
31   !!             -   !  2015-11  (J. Chanut) free surface simplification (remove filtered free surface)
32   !!            4.0  !  2017-05  (G. Madec)  introduction of the vertical physics manager (zdfphy)
33   !!----------------------------------------------------------------------
34
35   !!----------------------------------------------------------------------
36   !!   stp             : OPA system time-stepping
37   !!----------------------------------------------------------------------
38   USE step_oce         ! time stepping definition modules
39   !
40   USE iom              ! xIOs server
41
42   IMPLICIT NONE
43   PRIVATE
44
45   PUBLIC   stp   ! called by nemogcm.F90
46   PUBLIC   update_pointers ! called by nemo_init
47
48   !!----------------------------------------------------------------------
49   !! time level indices
50   !!----------------------------------------------------------------------
51   INTEGER, PUBLIC :: Nbb, Nnn, Naa, Nrhs          !! used by nemo_init
52
53   !!----------------------------------------------------------------------
54   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
55   !! $Id$
56   !! Software governed by the CeCILL license (see ./LICENSE)
57   !!----------------------------------------------------------------------
58CONTAINS
59
60#if defined key_agrif
61   RECURSIVE SUBROUTINE stp( )
62      INTEGER             ::   kstp   ! ocean time-step index
63#else
64   SUBROUTINE stp( kstp )
65      INTEGER, INTENT(in) ::   kstp   ! ocean time-step index
66#endif
67      !!----------------------------------------------------------------------
68      !!                     ***  ROUTINE stp  ***
69      !!
70      !! ** Purpose : - Time stepping of OPA  (momentum and active tracer eqs.)
71      !!              - Time stepping of SI3 (dynamic and thermodynamic eqs.)
72      !!              - Time stepping of TRC  (passive tracer eqs.)
73      !!
74      !! ** Method  : -1- Update forcings and data
75      !!              -2- Update ocean physics
76      !!              -3- Compute the t and s trends
77      !!              -4- Update t and s
78      !!              -5- Compute the momentum trends
79      !!              -6- Update the horizontal velocity
80      !!              -7- Compute the diagnostics variables (rd,N2, hdiv,w)
81      !!              -8- Outputs and diagnostics
82      !!----------------------------------------------------------------------
83      INTEGER ::   ji, jj, jk   ! dummy loop indice
84      INTEGER ::   indic        ! error indicator if < 0
85!!gm kcall can be removed, I guess
86      INTEGER ::   kcall        ! optional integer argument (dom_vvl_sf_nxt)
87      !! ---------------------------------------------------------------------
88#if defined key_agrif
89      kstp = nit000 + Agrif_Nb_Step()
90      IF( lk_agrif_debug ) THEN
91         IF( Agrif_Root() .and. lwp)   WRITE(*,*) '---'
92         IF(lwp)   WRITE(*,*) 'Grid Number', Agrif_Fixed(),' time step ', kstp, 'int tstep', Agrif_NbStepint()
93      ENDIF
94      IF( kstp == nit000 + 1 )   lk_agrif_fstep = .FALSE.
95# if defined key_iomput
96      IF( Agrif_Nbstepint() == 0 )   CALL iom_swap( cxios_context )
97# endif
98#endif
99      !
100      IF( ln_timing )   CALL timing_start('stp')
101      !
102      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
103      ! update I/O and calendar
104      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
105                             indic = 0                ! reset to no error condition
106                             
107      IF( kstp == nit000 ) THEN                       ! initialize IOM context (must be done after nemo_init for AGRIF+XIOS+OASIS)
108                             CALL iom_init(      cxios_context          )  ! for model grid (including passible AGRIF zoom)
109         IF( ln_crs      )   CALL iom_init( TRIM(cxios_context)//"_crs" )  ! for coarse grid
110      ENDIF
111      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
112                             CALL iom_setkt( kstp - nit000 + 1,      cxios_context          )   ! tell IOM we are at time step kstp
113      IF( ln_crs         )   CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" )   ! tell IOM we are at time step kstp
114
115      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
116      ! Update external forcing (tides, open boundaries, and surface boundary condition (including sea-ice)
117      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
118      IF( ln_tide    )   CALL sbc_tide( kstp )                   ! update tide potential
119      IF( ln_apr_dyn )   CALL sbc_apr ( kstp )                   ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib)
120      IF( ln_bdy     )   CALL bdy_dta ( kstp, time_offset=+1 )   ! update dynamic & tracer data at open boundaries
121                         CALL sbc     ( kstp, Nbb, Nnn )         ! Sea Boundary Condition (including sea-ice)
122
123      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
124      ! Update stochastic parameters and random T/S fluctuations
125      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
126      IF( ln_sto_eos ) CALL sto_par( kstp )          ! Stochastic parameters
127      IF( ln_sto_eos ) CALL sto_pts( tsn  )          ! Random T/S fluctuations
128
129      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
130      ! Ocean physics update
131      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
132      !  THERMODYNAMICS
133                         CALL eos_rab( tsb, rab_b )       ! before local thermal/haline expension ratio at T-points
134                         CALL eos_rab( tsn, rab_n )       ! now    local thermal/haline expension ratio at T-points
135                         CALL bn2    ( tsb, rab_b, rn2b ) ! before Brunt-Vaisala frequency
136                         CALL bn2    ( tsn, rab_n, rn2  ) ! now    Brunt-Vaisala frequency
137
138      !  VERTICAL PHYSICS
139                         CALL zdf_phy( kstp, Nbb, Nnn, Nrhs )   ! vertical physics update (top/bot drag, avt, avs, avm + MLD)
140
141      !  LATERAL  PHYSICS
142      !
143      IF( l_ldfslp ) THEN                             ! slope of lateral mixing
144                         CALL eos( tsb, rhd, gdept_0(:,:,:) )               ! before in situ density
145
146         IF( ln_zps .AND. .NOT. ln_isfcav)                               &
147            &            CALL zps_hde    ( kstp, jpts, tsb, gtsu, gtsv,  &  ! Partial steps: before horizontal gradient
148            &                                          rhd, gru , grv    )  ! of t, s, rd at the last ocean level
149
150         IF( ln_zps .AND.       ln_isfcav)                               &
151            &            CALL zps_hde_isf( kstp, jpts, tsb, gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF)
152            &                                          rhd, gru , grv , grui, grvi   )  ! of t, s, rd at the first ocean level
153         IF( ln_traldf_triad ) THEN
154                         CALL ldf_slp_triad( kstp, Nbb, Nnn )             ! before slope for triad operator
155         ELSE     
156                         CALL ldf_slp     ( kstp, rhd, rn2b, Nbb, Nnn )   ! before slope for standard operator
157         ENDIF
158      ENDIF
159      !                                                                        ! eddy diffusivity coeff.
160      IF( l_ldftra_time .OR. l_ldfeiv_time )   CALL ldf_tra( kstp, Nbb, Nnn )  !       and/or eiv coeff.
161      IF( l_ldfdyn_time                    )   CALL ldf_dyn( kstp, Nbb )       ! eddy viscosity coeff.
162
163      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
164      !  Ocean dynamics : hdiv, ssh, e3, u, v, w
165      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
166
167                            CALL ssh_nxt       ( kstp, Nnn )  ! after ssh (includes call to div_hor)
168      IF( .NOT.ln_linssh )  CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors
169                            CALL wzv           ( kstp )  ! now cross-level velocity
170      IF( ln_zad_Aimp )     CALL wAimp         ( kstp )  ! Adaptive-implicit vertical advection partitioning
171                            CALL eos    ( tsn, rhd, rhop, gdept_n(:,:,:) )  ! now in situ density for hpg computation
172                           
173!!jc: fs simplification
174!!jc: lines below are useless if ln_linssh=F. Keep them here (which maintains a bug if ln_linssh=T and ln_zps=T, cf ticket #1636)
175!!                                         but ensures reproductible results
176!!                                         with previous versions using split-explicit free surface         
177            IF( ln_zps .AND. .NOT. ln_isfcav )                               &
178               &            CALL zps_hde    ( kstp, jpts, tsn, gtsu, gtsv,   &  ! Partial steps: before horizontal gradient
179               &                                          rhd, gru , grv     )  ! of t, s, rd at the last ocean level
180            IF( ln_zps .AND.       ln_isfcav )                                          &
181               &            CALL zps_hde_isf( kstp, jpts, tsn, gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF)
182               &                                          rhd, gru , grv , grui, grvi   )  ! of t, s, rd at the first ocean level
183!!jc: fs simplification
184                           
185                         ua(:,:,:) = 0._wp            ! set dynamics trends to zero
186                         va(:,:,:) = 0._wp
187
188      IF(  lk_asminc .AND. ln_asmiau .AND. ln_dyninc )   &
189               &         CALL dyn_asm_inc   ( kstp )  ! apply dynamics assimilation increment
190      IF( ln_bdy     )   CALL bdy_dyn3d_dmp ( kstp )  ! bdy damping trends
191#if defined key_agrif
192      IF(.NOT. Agrif_Root())  & 
193               &         CALL Agrif_Sponge_dyn        ! momentum sponge
194#endif
195                         CALL dyn_adv( kstp, Nbb, Nnn      , uu, vv, Nrhs )  ! advection (VF or FF)   ==> RHS
196                         CALL dyn_vor( kstp,      Nnn      , uu, vv, Nrhs )  ! vorticity              ==> RHS
197                         CALL dyn_ldf( kstp, Nbb, Nnn      , uu, vv, Nrhs )  ! lateral mixing
198      IF( ln_zdfosm  )   CALL dyn_osm( kstp,      Nnn      , uu, vv, Nrhs )  ! OSMOSIS non-local velocity fluxes ==> RHS
199                         CALL dyn_hpg( kstp,      Nnn      , uu, vv, Nrhs )  ! horizontal gradient of Hydrostatic pressure
200                         CALL dyn_spg( kstp, Nbb, Nnn, Nrhs, uu, vv, ssh, uu_b, vv_b, Naa )  ! surface pressure gradient
201
202                                                      ! With split-explicit free surface, since now transports have been updated and ssha as well
203      IF( ln_dynspg_ts ) THEN                         ! vertical scale factors and vertical velocity need to be updated
204                            CALL div_hor    ( kstp, Nnn )         ! Horizontal divergence  (2nd call in time-split case)
205         IF(.NOT.ln_linssh) CALL dom_vvl_sf_nxt( kstp, kcall=2 )  ! after vertical scale factors (update depth average component)
206                            CALL wzv        ( kstp )              ! now cross-level velocity
207         IF( ln_zad_Aimp )  CALL wAimp      ( kstp )  ! Adaptive-implicit vertical advection partitioning
208      ENDIF
209     
210                         CALL dyn_zdf( kstp, Nbb, Nnn, Nrhs, uu, vv, Naa  )  ! vertical diffusion     ==> after
211
212      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
213      ! cool skin
214      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<     
215      IF ( ln_diurnal )  CALL diurnal_layers( kstp )
216     
217      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
218      ! diagnostics and outputs
219      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
220      IF( lk_floats  )   CALL flo_stp ( kstp )        ! drifting Floats
221      IF( ln_diacfl  )   CALL dia_cfl ( kstp )        ! Courant number diagnostics
222      IF( lk_diahth  )   CALL dia_hth ( kstp )        ! Thermocline depth (20 degres isotherm depth)
223      IF( lk_diadct  )   CALL dia_dct ( kstp )        ! Transports
224                         CALL dia_ar5 ( kstp )        ! ar5 diag
225      IF( lk_diaharm )   CALL dia_harm( kstp )        ! Tidal harmonic analysis
226                         CALL dia_wri ( kstp )        ! ocean model: outputs
227      !
228      IF( ln_crs     )   CALL crs_fld       ( kstp )  ! ocean model: online field coarsening & output
229     
230#if defined key_top
231      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
232      ! Passive Tracer Model
233      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
234                         CALL trc_stp       ( kstp, Nbb, Nnn, Nrhs, Naa )  ! time-stepping
235#endif
236
237      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
238      ! Active tracers                             
239      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
240                         tsa(:,:,:,:) = 0._wp         ! set tracer trends to zero
241
242      IF(  lk_asminc .AND. ln_asmiau .AND. &
243         & ln_trainc )   CALL tra_asm_inc   ( kstp )  ! apply tracer assimilation increment
244                         CALL tra_sbc       ( kstp, Nnn, Nrhs )  ! surface boundary condition
245      IF( ln_traqsr  )   CALL tra_qsr       ( kstp, Nnn, Nrhs )  ! penetrative solar radiation qsr
246      IF( ln_trabbc  )   CALL tra_bbc       ( kstp, Nnn, Nrhs )  ! bottom heat flux
247      IF( ln_trabbl  )   CALL tra_bbl       ( kstp, Nnn, Nrhs )  ! advective (and/or diffusive) bottom boundary layer scheme
248      IF( ln_tradmp  )   CALL tra_dmp       ( kstp, Nnn, Nrhs )  ! internal damping trends
249      IF( ln_bdy     )   CALL bdy_tra_dmp   ( kstp )  ! bdy damping trends
250#if defined key_agrif
251      IF(.NOT. Agrif_Root())  & 
252               &         CALL Agrif_Sponge_tra        ! tracers sponge
253#endif
254                         CALL tra_adv( kstp, Nbb, Nnn      , ts, Nrhs )  ! hor. + vert. advection  ==> RHS
255      IF( ln_zdfosm  )   CALL tra_osm( kstp,      Nnn      , ts, Nrhs )  ! OSMOSIS non-local tracer fluxes ==> RHS
256      IF( lrst_oce .AND. ln_zdfosm ) &
257           &             CALL osm_rst( kstp, Nnn, 'WRITE' )! write OSMOSIS outputs + wn (so must do here) to restarts
258                         CALL tra_ldf( kstp, Nnn, Nrhs )   ! lateral mixing
259
260!!gm : why CALL to dia_ptr has been moved here??? (use trends info?)
261      IF( ln_diaptr  )   CALL dia_ptr                 ! Poleward adv/ldf TRansports diagnostics
262!!gm
263                         CALL tra_zdf( kstp, Nbb, Nnn, Nrhs, ts, Naa  )  ! vert. mixing & after tracer   ==> after
264      IF( ln_zdfnpc  )   CALL tra_npc( kstp,      Nnn, Nrhs           )  ! update after fields by non-penetrative convection
265
266      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
267      ! Set boundary conditions and Swap
268      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
269!!jc1: For agrif, it would be much better to finalize tracers/momentum here (e.g. bdy conditions) and move the swap
270!!    (and time filtering) after Agrif update. Then restart would be done after and would contain updated fields.
271!!    If so:
272!!    (i) no need to call agrif update at initialization time
273!!    (ii) no need to update "before" fields
274!!
275!!    Apart from creating new tra_swp/dyn_swp routines, this however:
276!!    (i) makes boundary conditions at initialization time computed from updated fields which is not the case between
277!!    two restarts => restartability issue. One can circumvent this, maybe, by assuming "interface separation",
278!!    e.g. a shift of the feedback interface inside child domain.
279!!    (ii) requires that all restart outputs of updated variables by agrif (e.g. passive tracers/tke/barotropic arrays) are done at the same
280!!    place.
281!!
282!!jc2: dynnxt must be the latest call. e3t_b are indeed updated in that routine
283                         CALL tra_nxt       ( kstp, Nnn, Nrhs )  ! finalize (bcs) tracer fields at next time step and swap
284                         CALL dyn_nxt       ( kstp, Nnn )        ! finalize (bcs) velocities at next time step and swap (always called after tra_nxt)
285                         CALL ssh_swp       ( kstp )  ! swap of sea surface height
286      IF(.NOT.ln_linssh) CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors
287      !
288      IF( ln_diahsb  )   CALL dia_hsb       ( kstp )  ! - ML - global conservation diagnostics
289
290!!gm : This does not only concern the dynamics ==>>> add a new title
291!!gm2: why ouput restart before AGRIF update?
292!!
293!!jc: That would be better, but see comment above
294!!
295      IF( lrst_oce   )   CALL rst_write    ( kstp, Nbb, Nnn )   ! write output ocean restart file
296      IF( ln_sto_eos )   CALL sto_rst_write( kstp )   ! write restart file for stochastic parameters
297
298#if defined key_agrif
299      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
300      ! AGRIF
301      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<     
302                         CALL Agrif_Integrate_ChildGrids( stp )  ! allows to finish all the Child Grids before updating
303
304                         IF( Agrif_NbStepint() == 0 ) CALL Agrif_update_all( ) ! Update all components
305#endif
306      IF( ln_diaobs  )   CALL dia_obs      ( kstp, Nnn )      ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
307
308      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
309      ! Control
310      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
311                         CALL stp_ctl      ( kstp, indic )
312                         
313      IF( kstp == nit000 ) THEN                          ! 1st time step only
314                                        CALL iom_close( numror )   ! close input  ocean restart file
315         IF(lwm)                        CALL FLUSH    ( numond )   ! flush output namelist oce
316         IF(lwm .AND. numoni /= -1 )    CALL FLUSH    ( numoni )   ! flush output namelist ice (if exist)
317      ENDIF
318
319      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
320      ! Coupled mode
321      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
322!!gm why lk_oasis and not lk_cpl ????
323      IF( lk_oasis   )   CALL sbc_cpl_snd( kstp, Nnn )     ! coupled mode : field exchanges
324      !
325#if defined key_iomput
326      IF( kstp == nitend .OR. indic < 0 ) THEN
327                      CALL iom_context_finalize(      cxios_context          ) ! needed for XIOS+AGRIF
328                      IF(lrxios) CALL iom_context_finalize(      crxios_context          )
329         IF( ln_crs ) CALL iom_context_finalize( trim(cxios_context)//"_crs" ) !
330      ENDIF
331#endif
332      !
333      IF( ln_timing )   CALL timing_stop('stp')
334      !
335   END SUBROUTINE stp
336   
337   SUBROUTINE update_pointers( Kbb, Kmm, Kaa )
338      !!----------------------------------------------------------------------
339      !!                     ***  ROUTINE update_pointers  ***
340      !!
341      !! ** Purpose :   Associate temporary pointer arrays.
342      !!                For IMMERSE development phase only - to be deleted
343      !!
344      !! ** Method  :
345      !!----------------------------------------------------------------------
346      INTEGER, INTENT( in ) :: Kbb, Kmm, Kaa ! time level indices
347
348      ub => uu(:,:,:,Kbb); un => uu(:,:,:,Kmm); ua => uu(:,:,:,Kaa)
349      vb => vv(:,:,:,Kbb); vn => vv(:,:,:,Kmm); va => vv(:,:,:,Kaa)
350      wn => ww(:,:,:)
351      hdivn => hdiv(:,:,:)
352
353      sshb =>  ssh(:,:,Kbb); sshn =>  ssh(:,:,Kmm); ssha =>  ssh(:,:,Kaa)
354      ub_b => uu_b(:,:,Kbb); un_b => uu_b(:,:,Kmm); ua_b => uu_b(:,:,Kaa)
355      vb_b => vv_b(:,:,Kbb); vn_b => vv_b(:,:,Kmm); va_b => vv_b(:,:,Kaa)
356
357      tsb => ts(:,:,:,:,Kbb); tsn => ts(:,:,:,:,Kmm); tsa => ts(:,:,:,:,Kaa)
358
359      e3t_b => e3t(:,:,:,Kbb); e3t_n => e3t(:,:,:,Kmm); e3t_a => e3t(:,:,:,Kaa)
360      e3u_b => e3u(:,:,:,Kbb); e3u_n => e3u(:,:,:,Kmm); e3u_a => e3u(:,:,:,Kaa)
361      e3v_b => e3v(:,:,:,Kbb); e3v_n => e3v(:,:,:,Kmm); e3v_a => e3v(:,:,:,Kaa)
362
363      e3f_n => e3f(:,:,:)
364
365      e3w_b  => e3w (:,:,:,Kbb); e3w_n  => e3w (:,:,:,Kmm)
366      e3uw_b => e3uw(:,:,:,Kbb); e3uw_n => e3uw(:,:,:,Kmm)
367      e3vw_b => e3vw(:,:,:,Kbb); e3vw_n => e3vw(:,:,:,Kmm)
368
369      gdept_b => gdept(:,:,:,Kbb); gdept_n => gdept(:,:,:,Kmm) 
370      gdepw_b => gdepw(:,:,:,Kbb); gdepw_n => gdepw(:,:,:,Kmm) 
371      gde3w_n => gde3w(:,:,:)
372
373   END SUBROUTINE update_pointers
374
375   !!======================================================================
376END MODULE step
Note: See TracBrowser for help on using the repository browser.