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.
dtadyn.F90 in NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps_rewrite_time_filterswap/src/OFF – NEMO

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps_rewrite_time_filterswap/src/OFF/dtadyn.F90 @ 11427

Last change on this file since 11427 was 11427, checked in by davestorkey, 5 years ago

dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps_rewrite_time_filterswap :
Restore independent time level indices for TOP in order to enable the option for
TOP to have a different timestep to OCE (nn_dttrc > 1). But note that this version of
the code only works for nn_dttrc=1. Also sort out the time-level swapping for OFF.
This commit passes the GYRE_PISCES and ORCA2_OFF_PISCES tests but fails restartability
and bit-comparison with the control for ORCA2_ICE_PISCES.

  • Property svn:keywords set to Id
File size: 42.3 KB
Line 
1MODULE dtadyn
2   !!======================================================================
3   !!                       ***  MODULE  dtadyn  ***
4   !! Off-line : interpolation of the physical fields
5   !!======================================================================
6   !! History :   OPA  ! 1992-01 (M. Imbard) Original code
7   !!             8.0  ! 1998-04 (L.Bopp MA Foujols) slopes for isopyc.
8   !!              -   ! 1998-05 (L. Bopp) read output of coupled run
9   !!             8.2  ! 2001-01 (M. Levy et M. Benjelloul) add netcdf FORMAT
10   !!   NEMO      1.0  ! 2005-03 (O. Aumont and A. El Moussaoui) F90
11   !!              -   ! 2005-12 (C. Ethe) Adapted for DEGINT
12   !!             3.0  ! 2007-06 (C. Ethe) use of iom module
13   !!             3.3  ! 2010-11 (C. Ethe) Full reorganization of the off-line: phasing with the on-line
14   !!             3.4  ! 2011-05 (C. Ethe) Use of fldread
15   !!----------------------------------------------------------------------
16
17   !!----------------------------------------------------------------------
18   !!   dta_dyn_init : initialization, namelist read, and SAVEs control
19   !!   dta_dyn      : Interpolation of the fields
20   !!----------------------------------------------------------------------
21   USE oce             ! ocean dynamics and tracers variables
22   USE c1d             ! 1D configuration: lk_c1d
23   USE dom_oce         ! ocean domain: variables
24   USE domvvl          ! variable volume
25   USE zdf_oce         ! ocean vertical physics: variables
26   USE sbc_oce         ! surface module: variables
27   USE trc_oce         ! share ocean/biogeo variables
28   USE phycst          ! physical constants
29   USE trabbl          ! active tracer: bottom boundary layer
30   USE ldfslp          ! lateral diffusion: iso-neutral slopes
31   USE sbcrnf          ! river runoffs
32   USE ldftra          ! ocean tracer   lateral physics
33   USE zdfmxl          ! vertical physics: mixed layer depth
34   USE eosbn2          ! equation of state - Brunt Vaisala frequency
35   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
36   USE zpshde          ! z-coord. with partial steps: horizontal derivatives
37   USE in_out_manager  ! I/O manager
38   USE iom             ! I/O library
39   USE lib_mpp         ! distributed memory computing library
40   USE prtctl          ! print control
41   USE fldread         ! read input fields
42   USE timing          ! Timing
43   USE trc, ONLY : ln_rsttr, numrtr, numrtw, lrst_trc
44
45   IMPLICIT NONE
46   PRIVATE
47
48   PUBLIC   dta_dyn_init       ! called by nemo_init
49   PUBLIC   dta_dyn            ! called by nemo_gcm
50   PUBLIC   dta_dyn_sed_init   ! called by nemo_init
51   PUBLIC   dta_dyn_sed        ! called by nemo_gcm
52   PUBLIC   dta_dyn_atf        ! called by nemo_gcm
53   PUBLIC   dta_dyn_sf_interp  ! called by nemo_gcm
54
55   CHARACTER(len=100) ::   cn_dir          !: Root directory for location of ssr files
56   LOGICAL            ::   ln_dynrnf       !: read runoff data in file (T) or set to zero (F)
57   LOGICAL            ::   ln_dynrnf_depth       !: read runoff data in file (T) or set to zero (F)
58   REAL(wp)           ::   fwbcorr
59
60
61   INTEGER  , PARAMETER ::   jpfld = 20     ! maximum number of fields to read
62   INTEGER  , SAVE      ::   jf_tem         ! index of temperature
63   INTEGER  , SAVE      ::   jf_sal         ! index of salinity
64   INTEGER  , SAVE      ::   jf_uwd         ! index of u-transport
65   INTEGER  , SAVE      ::   jf_vwd         ! index of v-transport
66   INTEGER  , SAVE      ::   jf_wwd         ! index of v-transport
67   INTEGER  , SAVE      ::   jf_avt         ! index of Kz
68   INTEGER  , SAVE      ::   jf_mld         ! index of mixed layer deptht
69   INTEGER  , SAVE      ::   jf_emp         ! index of water flux
70   INTEGER  , SAVE      ::   jf_empb        ! index of water flux
71   INTEGER  , SAVE      ::   jf_qsr         ! index of solar radiation
72   INTEGER  , SAVE      ::   jf_wnd         ! index of wind speed
73   INTEGER  , SAVE      ::   jf_ice         ! index of sea ice cover
74   INTEGER  , SAVE      ::   jf_rnf         ! index of river runoff
75   INTEGER  , SAVE      ::   jf_fmf         ! index of downward salt flux
76   INTEGER  , SAVE      ::   jf_ubl         ! index of u-bbl coef
77   INTEGER  , SAVE      ::   jf_vbl         ! index of v-bbl coef
78   INTEGER  , SAVE      ::   jf_div         ! index of e3t
79
80
81   TYPE(FLD), ALLOCATABLE, SAVE, DIMENSION(:) :: sf_dyn  ! structure of input fields (file informations, fields read)
82   !                                               !
83   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: uslpdta    ! zonal isopycnal slopes
84   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: vslpdta    ! meridional isopycnal slopes
85   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpidta   ! zonal diapycnal slopes
86   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpjdta   ! meridional diapycnal slopes
87
88   INTEGER, SAVE  :: nprevrec, nsecdyn
89
90   !!----------------------------------------------------------------------
91   !! NEMO/OFF 4.0 , NEMO Consortium (2018)
92   !! $Id$
93   !! Software governed by the CeCILL license (see ./LICENSE)
94   !!----------------------------------------------------------------------
95CONTAINS
96
97   SUBROUTINE dta_dyn( kt, Kbb, Kmm, Kaa )
98      !!----------------------------------------------------------------------
99      !!                  ***  ROUTINE dta_dyn  ***
100      !!
101      !! ** Purpose :  Prepares dynamics and physics fields from a NEMO run
102      !!               for an off-line simulation of passive tracers
103      !!
104      !! ** Method : calculates the position of data
105      !!             - computes slopes if needed
106      !!             - interpolates data if needed
107      !!----------------------------------------------------------------------
108      INTEGER, INTENT(in) ::   kt             ! ocean time-step index
109      INTEGER, INTENT(in) ::   Kbb, Kmm, Kaa  ! ocean time level indices
110      !
111      INTEGER             ::   ji, jj, jk
112      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   zemp
113      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   zhdivtr
114      !!----------------------------------------------------------------------
115      !
116      IF( ln_timing )   CALL timing_start( 'dta_dyn')
117      !
118      nsecdyn = nsec_year + nsec1jan000   ! number of seconds between Jan. 1st 00h of nit000 year and the middle of time step
119      !
120      IF( kt == nit000 ) THEN    ;    nprevrec = 0
121      ELSE                       ;    nprevrec = sf_dyn(jf_tem)%nrec_a(2)
122      ENDIF
123      CALL fld_read( kt, 1, sf_dyn )      !=  read data at kt time step   ==!
124      !
125      IF( l_ldfslp .AND. .NOT.lk_c1d )   CALL  dta_dyn_slp( kt, Kbb, Kmm )    ! Computation of slopes
126      !
127      ts(:,:,:,jp_tem,Kmm) = sf_dyn(jf_tem)%fnow(:,:,:)  * tmask(:,:,:)    ! temperature
128      ts(:,:,:,jp_sal,Kmm) = sf_dyn(jf_sal)%fnow(:,:,:)  * tmask(:,:,:)    ! salinity
129      wndm(:,:)         = sf_dyn(jf_wnd)%fnow(:,:,1)  * tmask(:,:,1)    ! wind speed - needed for gas exchange
130      fmmflx(:,:)       = sf_dyn(jf_fmf)%fnow(:,:,1)  * tmask(:,:,1)    ! downward salt flux (v3.5+)
131      fr_i(:,:)         = sf_dyn(jf_ice)%fnow(:,:,1)  * tmask(:,:,1)    ! Sea-ice fraction
132      qsr (:,:)         = sf_dyn(jf_qsr)%fnow(:,:,1)  * tmask(:,:,1)    ! solar radiation
133      emp (:,:)         = sf_dyn(jf_emp)%fnow(:,:,1)  * tmask(:,:,1)    ! E-P
134      IF( ln_dynrnf ) THEN
135         rnf (:,:)      = sf_dyn(jf_rnf)%fnow(:,:,1) * tmask(:,:,1)    ! E-P
136         IF( ln_dynrnf_depth .AND. .NOT. ln_linssh )    CALL  dta_dyn_hrnf(Kmm)
137      ENDIF
138      !
139      uu(:,:,:,Kmm)        = sf_dyn(jf_uwd)%fnow(:,:,:) * umask(:,:,:)    ! effective u-transport
140      vv(:,:,:,Kmm)        = sf_dyn(jf_vwd)%fnow(:,:,:) * vmask(:,:,:)    ! effective v-transport
141      ww(:,:,:)        = sf_dyn(jf_wwd)%fnow(:,:,:) * tmask(:,:,:)    ! effective v-transport
142      !
143      IF( .NOT.ln_linssh ) THEN
144         ALLOCATE( zemp(jpi,jpj) , zhdivtr(jpi,jpj,jpk) )
145         zhdivtr(:,:,:) = sf_dyn(jf_div)%fnow(:,:,:)  * tmask(:,:,:)    ! effective u-transport
146         emp_b  (:,:)   = sf_dyn(jf_empb)%fnow(:,:,1) * tmask(:,:,1)    ! E-P
147         zemp   (:,:)   = ( 0.5_wp * ( emp(:,:) + emp_b(:,:) ) + rnf(:,:) + fwbcorr ) * tmask(:,:,1)
148         CALL dta_dyn_ssh( kt, zhdivtr, ssh(:,:,Kbb), zemp, ssh(:,:,Kaa), e3t(:,:,:,Kaa) )  !=  ssh, vertical scale factor & vertical transport
149         DEALLOCATE( zemp , zhdivtr )
150         !                                           Write in the tracer restart file
151         !                                          *********************************
152         IF( lrst_trc ) THEN
153            IF(lwp) WRITE(numout,*)
154            IF(lwp) WRITE(numout,*) 'dta_dyn_ssh : ssh field written in tracer restart file at it= ', kt,' date= ', ndastp
155            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
156            CALL iom_rstput( kt, nitrst, numrtw, 'sshn', ssh(:,:,Kaa) )
157            CALL iom_rstput( kt, nitrst, numrtw, 'sshb', ssh(:,:,Kmm) )
158         ENDIF
159      ENDIF
160      !
161      CALL eos    ( ts(:,:,:,:,Kmm), rhd, rhop, gdept_0(:,:,:) ) ! In any case, we need rhop
162      CALL eos_rab( ts(:,:,:,:,Kmm), rab_n, Kmm )       ! now    local thermal/haline expension ratio at T-points
163      CALL bn2    ( ts(:,:,:,:,Kmm), rab_n, rn2, Kmm )  ! before Brunt-Vaisala frequency need for zdfmxl
164
165      rn2b(:,:,:) = rn2(:,:,:)         ! needed for zdfmxl
166      CALL zdf_mxl( kt, Kmm )          ! In any case, we need mxl
167      !
168      hmld(:,:)       = sf_dyn(jf_mld)%fnow(:,:,1) * tmask(:,:,1)    ! mixed layer depht
169      avt(:,:,:)      = sf_dyn(jf_avt)%fnow(:,:,:) * tmask(:,:,:)    ! vertical diffusive coefficient
170      avs(:,:,:)      = avt(:,:,:)
171      !
172      IF( ln_trabbl .AND. .NOT.lk_c1d ) THEN       ! diffusive Bottom boundary layer param
173         ahu_bbl(:,:) = sf_dyn(jf_ubl)%fnow(:,:,1) * umask(:,:,1)    ! bbl diffusive coef
174         ahv_bbl(:,:) = sf_dyn(jf_vbl)%fnow(:,:,1) * vmask(:,:,1)
175      ENDIF
176      !
177      !
178      CALL eos( ts(:,:,:,:,Kmm), rhd, rhop, gdept_0(:,:,:) ) ! In any case, we need rhop
179      !
180      IF(ln_ctl) THEN                  ! print control
181         CALL prt_ctl(tab3d_1=ts(:,:,:,jp_tem,Kmm), clinfo1=' tn      - : ', mask1=tmask,  kdim=jpk   )
182         CALL prt_ctl(tab3d_1=ts(:,:,:,jp_sal,Kmm), clinfo1=' sn      - : ', mask1=tmask,  kdim=jpk   )
183         CALL prt_ctl(tab3d_1=uu(:,:,:,Kmm)               , clinfo1=' uu(:,:,:,Kmm)      - : ', mask1=umask,  kdim=jpk   )
184         CALL prt_ctl(tab3d_1=vv(:,:,:,Kmm)               , clinfo1=' vv(:,:,:,Kmm)      - : ', mask1=vmask,  kdim=jpk   )
185         CALL prt_ctl(tab3d_1=ww               , clinfo1=' ww      - : ', mask1=tmask,  kdim=jpk   )
186         CALL prt_ctl(tab3d_1=avt              , clinfo1=' kz      - : ', mask1=tmask,  kdim=jpk   )
187         CALL prt_ctl(tab3d_1=uslp             , clinfo1=' slp  - u : ', tab3d_2=vslp, clinfo2=' v : ', kdim=jpk)
188         CALL prt_ctl(tab3d_1=wslpi            , clinfo1=' slp  - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk)
189      ENDIF
190      !
191      IF( ln_timing )   CALL timing_stop( 'dta_dyn')
192      !
193   END SUBROUTINE dta_dyn
194
195
196   SUBROUTINE dta_dyn_init( Kbb, Kmm, Kaa )
197      !!----------------------------------------------------------------------
198      !!                  ***  ROUTINE dta_dyn_init  ***
199      !!
200      !! ** Purpose :   Initialisation of the dynamical data     
201      !! ** Method  : - read the data namdta_dyn namelist
202      !!----------------------------------------------------------------------
203      INTEGER, INTENT(in) ::   Kbb, Kmm, Kaa  ! ocean time level indices
204      !
205      INTEGER  :: ierr, ierr0, ierr1, ierr2, ierr3   ! return error code
206      INTEGER  :: ifpr                               ! dummy loop indice
207      INTEGER  :: jfld                               ! dummy loop arguments
208      INTEGER  :: inum, idv, idimv                   ! local integer
209      INTEGER  :: ios                                ! Local integer output status for namelist read
210      INTEGER  :: ji, jj, jk
211      REAL(wp) :: zcoef
212      INTEGER  :: nkrnf_max
213      REAL(wp) :: hrnf_max
214      !!
215      CHARACTER(len=100)            ::  cn_dir        !   Root directory for location of core files
216      TYPE(FLD_N), DIMENSION(jpfld) ::  slf_d         ! array of namelist informations on the fields to read
217      TYPE(FLD_N) :: sn_uwd, sn_vwd, sn_wwd, sn_empb, sn_emp  ! informations about the fields to be read
218      TYPE(FLD_N) :: sn_tem , sn_sal , sn_avt   !   "                 "
219      TYPE(FLD_N) :: sn_mld, sn_qsr, sn_wnd , sn_ice , sn_fmf   !   "               "
220      TYPE(FLD_N) :: sn_ubl, sn_vbl, sn_rnf    !   "              "
221      TYPE(FLD_N) :: sn_div  ! informations about the fields to be read
222      !!
223      NAMELIST/namdta_dyn/cn_dir, ln_dynrnf, ln_dynrnf_depth,  fwbcorr, &
224         &                sn_uwd, sn_vwd, sn_wwd, sn_emp,               &
225         &                sn_avt, sn_tem, sn_sal, sn_mld , sn_qsr ,     &
226         &                sn_wnd, sn_ice, sn_fmf,                       &
227         &                sn_ubl, sn_vbl, sn_rnf,                       &
228         &                sn_empb, sn_div 
229      !!----------------------------------------------------------------------
230      !
231      REWIND( numnam_ref )              ! Namelist namdta_dyn in reference namelist : Offline: init. of dynamical data
232      READ  ( numnam_ref, namdta_dyn, IOSTAT = ios, ERR = 901)
233901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namdta_dyn in reference namelist', lwp )
234      REWIND( numnam_cfg )              ! Namelist namdta_dyn in configuration namelist : Offline: init. of dynamical data
235      READ  ( numnam_cfg, namdta_dyn, IOSTAT = ios, ERR = 902 )
236902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namdta_dyn in configuration namelist', lwp )
237      IF(lwm) WRITE ( numond, namdta_dyn )
238      !                                         ! store namelist information in an array
239      !                                         ! Control print
240      IF(lwp) THEN
241         WRITE(numout,*)
242         WRITE(numout,*) 'dta_dyn : offline dynamics '
243         WRITE(numout,*) '~~~~~~~ '
244         WRITE(numout,*) '   Namelist namdta_dyn'
245         WRITE(numout,*) '      runoffs option enabled (T) or not (F)            ln_dynrnf        = ', ln_dynrnf
246         WRITE(numout,*) '      runoffs is spread in vertical                    ln_dynrnf_depth  = ', ln_dynrnf_depth
247         WRITE(numout,*) '      annual global mean of empmr for ssh correction   fwbcorr          = ', fwbcorr
248         WRITE(numout,*)
249      ENDIF
250      !
251      jf_uwd  = 1     ;   jf_vwd  = 2    ;   jf_wwd = 3    ;   jf_emp = 4    ;   jf_avt = 5
252      jf_tem  = 6     ;   jf_sal  = 7    ;   jf_mld = 8    ;   jf_qsr = 9
253      jf_wnd  = 10    ;   jf_ice  = 11   ;   jf_fmf = 12   ;   jfld   = jf_fmf
254      !
255      slf_d(jf_uwd)  = sn_uwd    ;   slf_d(jf_vwd)  = sn_vwd   ;   slf_d(jf_wwd) = sn_wwd
256      slf_d(jf_emp)  = sn_emp    ;   slf_d(jf_avt)  = sn_avt
257      slf_d(jf_tem)  = sn_tem    ;   slf_d(jf_sal)  = sn_sal   ;   slf_d(jf_mld) = sn_mld
258      slf_d(jf_qsr)  = sn_qsr    ;   slf_d(jf_wnd)  = sn_wnd   ;   slf_d(jf_ice) = sn_ice
259      slf_d(jf_fmf)  = sn_fmf
260      !
261      IF( .NOT.ln_linssh ) THEN
262               jf_div  = jfld + 1   ;         jf_empb  = jfld + 2    ;   jfld = jf_empb
263         slf_d(jf_div) = sn_div     ;   slf_d(jf_empb) = sn_empb
264      ENDIF
265      !
266      IF( ln_trabbl ) THEN
267               jf_ubl  = jfld + 1   ;         jf_vbl  = jfld + 2     ;   jfld = jf_vbl
268         slf_d(jf_ubl) = sn_ubl     ;   slf_d(jf_vbl) = sn_vbl
269      ENDIF
270      !
271      IF( ln_dynrnf ) THEN
272               jf_rnf  = jfld + 1   ;     jfld  = jf_rnf
273         slf_d(jf_rnf) = sn_rnf   
274      ELSE
275         rnf(:,:) = 0._wp
276      ENDIF
277
278      ALLOCATE( sf_dyn(jfld), STAT=ierr )         ! set sf structure
279      IF( ierr > 0 )  THEN
280         CALL ctl_stop( 'dta_dyn: unable to allocate sf structure' )   ;   RETURN
281      ENDIF
282      !                                         ! fill sf with slf_i and control print
283      CALL fld_fill( sf_dyn, slf_d, cn_dir, 'dta_dyn_init', 'Data in file', 'namdta_dyn' )
284      !
285      ! Open file for each variable to get his number of dimension
286      DO ifpr = 1, jfld
287         CALL fld_clopn( sf_dyn(ifpr), nyear, nmonth, nday )
288         idv   = iom_varid( sf_dyn(ifpr)%num , slf_d(ifpr)%clvar )        ! id of the variable sdjf%clvar
289         idimv = iom_file ( sf_dyn(ifpr)%num )%ndims(idv)                 ! number of dimension for variable sdjf%clvar
290         IF( sf_dyn(ifpr)%num /= 0 )   CALL iom_close( sf_dyn(ifpr)%num ) ! close file if already open
291         ierr1=0
292         IF( idimv == 3 ) THEN    ! 2D variable
293                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,1)    , STAT=ierr0 )
294            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,1,2)  , STAT=ierr1 )
295         ELSE                     ! 3D variable
296                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,jpk)  , STAT=ierr0 )
297            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,jpk,2), STAT=ierr1 )
298         ENDIF
299         IF( ierr0 + ierr1 > 0 ) THEN
300            CALL ctl_stop( 'dta_dyn_init : unable to allocate sf_dyn array structure' )   ;   RETURN
301         ENDIF
302      END DO
303      !
304      IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN                  ! slopes
305         IF( sf_dyn(jf_tem)%ln_tint ) THEN      ! time interpolation
306            ALLOCATE( uslpdta (jpi,jpj,jpk,2), vslpdta (jpi,jpj,jpk,2),    &
307            &         wslpidta(jpi,jpj,jpk,2), wslpjdta(jpi,jpj,jpk,2), STAT=ierr2 )
308            !
309            IF( ierr2 > 0 )  THEN
310               CALL ctl_stop( 'dta_dyn_init : unable to allocate slope arrays' )   ;   RETURN
311            ENDIF
312         ENDIF
313      ENDIF
314      !
315      IF( .NOT.ln_linssh ) THEN
316        IF( .NOT. sf_dyn(jf_uwd)%ln_clim .AND. ln_rsttr .AND.    &                     ! Restart: read in restart file
317           iom_varid( numrtr, 'sshn', ldstop = .FALSE. ) > 0 ) THEN
318           IF(lwp) WRITE(numout,*) ' ssh(:,:,Kmm) forcing fields read in the restart file for initialisation'
319           CALL iom_get( numrtr, jpdom_autoglo, 'sshn', ssh(:,:,Kmm)   )
320           CALL iom_get( numrtr, jpdom_autoglo, 'sshb', ssh(:,:,Kbb)   )
321        ELSE
322           IF(lwp) WRITE(numout,*) ' ssh(:,:,Kmm) forcing fields read in the restart file for initialisation'
323           CALL iom_open( 'restart', inum )
324           CALL iom_get( inum, jpdom_autoglo, 'sshn', ssh(:,:,Kmm)   )
325           CALL iom_get( inum, jpdom_autoglo, 'sshb', ssh(:,:,Kbb)   )
326           CALL iom_close( inum )                                        ! close file
327        ENDIF
328        !
329        DO jk = 1, jpkm1
330           e3t(:,:,jk,Kmm) = e3t_0(:,:,jk) * ( 1._wp + ssh(:,:,Kmm) * tmask(:,:,1) / ( ht_0(:,:) + 1.0 - tmask(:,:,1) ) )
331        ENDDO
332        e3t(:,:,jpk,Kaa) = e3t_0(:,:,jpk)
333
334        ! Horizontal scale factor interpolations
335        ! --------------------------------------
336        CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3u(:,:,:,Kmm), 'U' )
337        CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3v(:,:,:,Kmm), 'V' )
338
339        ! Vertical scale factor interpolations
340        ! ------------------------------------
341        CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3w(:,:,:,Kmm), 'W' )
342 
343        e3t(:,:,:,Kbb)  = e3t(:,:,:,Kmm)
344        e3u(:,:,:,Kbb)  = e3u(:,:,:,Kmm)
345        e3v(:,:,:,Kbb)  = e3v(:,:,:,Kmm)
346
347        ! t- and w- points depth
348        ! ----------------------
349        gdept(:,:,1,Kmm) = 0.5_wp * e3w(:,:,1,Kmm)
350        gdepw(:,:,1,Kmm) = 0.0_wp
351
352        DO jk = 2, jpk
353           DO jj = 1,jpj
354              DO ji = 1,jpi
355                !    zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))   ! 0 everywhere
356                !    tmask = wmask, ie everywhere expect at jk = mikt
357                                                                   ! 1 for jk =
358                                                                   ! mikt
359                 zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
360                 gdepw(ji,jj,jk,Kmm) = gdepw(ji,jj,jk-1,Kmm) + e3t(ji,jj,jk-1,Kmm)
361                 gdept(ji,jj,jk,Kmm) =      zcoef  * ( gdepw(ji,jj,jk  ,Kmm) + 0.5 * e3w(ji,jj,jk,Kmm))  &
362                     &                + (1-zcoef) * ( gdept(ji,jj,jk-1,Kmm) + e3w(ji,jj,jk,Kmm))
363              END DO
364           END DO
365        END DO
366
367        gdept(:,:,:,Kbb) = gdept(:,:,:,Kmm)
368        gdepw(:,:,:,Kbb) = gdepw(:,:,:,Kmm)
369        !
370      ENDIF
371      !
372      IF( ln_dynrnf .AND. ln_dynrnf_depth ) THEN       ! read depht over which runoffs are distributed
373         IF(lwp) WRITE(numout,*) 
374         IF(lwp) WRITE(numout,*) ' read in the file depht over which runoffs are distributed'
375         CALL iom_open ( "runoffs", inum )                           ! open file
376         CALL iom_get  ( inum, jpdom_data, 'rodepth', h_rnf )   ! read the river mouth array
377         CALL iom_close( inum )                                        ! close file
378         !
379         nk_rnf(:,:) = 0                               ! set the number of level over which river runoffs are applied
380         DO jj = 1, jpj
381            DO ji = 1, jpi
382               IF( h_rnf(ji,jj) > 0._wp ) THEN
383                  jk = 2
384                  DO WHILE ( jk /= mbkt(ji,jj) .AND. gdept_0(ji,jj,jk) < h_rnf(ji,jj) ) ;  jk = jk + 1
385                  END DO
386                  nk_rnf(ji,jj) = jk
387               ELSEIF( h_rnf(ji,jj) == -1._wp   ) THEN   ;  nk_rnf(ji,jj) = 1
388               ELSEIF( h_rnf(ji,jj) == -999._wp ) THEN   ;  nk_rnf(ji,jj) = mbkt(ji,jj)
389               ELSE
390                  CALL ctl_stop( 'sbc_rnf_init: runoff depth not positive, and not -999 or -1, rnf value in file fort.999'  )
391                  WRITE(999,*) 'ji, jj, h_rnf(ji,jj) :', ji, jj, h_rnf(ji,jj)
392               ENDIF
393            END DO
394         END DO
395         DO jj = 1, jpj                                ! set the associated depth
396            DO ji = 1, jpi
397               h_rnf(ji,jj) = 0._wp
398               DO jk = 1, nk_rnf(ji,jj)
399                  h_rnf(ji,jj) = h_rnf(ji,jj) + e3t(ji,jj,jk,Kmm)
400               END DO
401            END DO
402         END DO
403      ELSE                                       ! runoffs applied at the surface
404         nk_rnf(:,:) = 1
405         h_rnf (:,:) = e3t(:,:,1,Kmm)
406      ENDIF
407      nkrnf_max = MAXVAL( nk_rnf(:,:) )
408      hrnf_max = MAXVAL( h_rnf(:,:) )
409      IF( lk_mpp )  THEN
410         CALL mpp_max( 'dtadyn', nkrnf_max )                 ! max over the  global domain
411         CALL mpp_max( 'dtadyn', hrnf_max )                 ! max over the  global domain
412      ENDIF
413      IF(lwp) WRITE(numout,*) ' '
414      IF(lwp) WRITE(numout,*) ' max depht of runoff : ', hrnf_max,'    max level  : ', nkrnf_max
415      IF(lwp) WRITE(numout,*) ' '
416      !
417      CALL dta_dyn( nit000, Kbb, Kmm, Kaa )
418      !
419   END SUBROUTINE dta_dyn_init
420
421   SUBROUTINE dta_dyn_sed( kt, Kmm )
422      !!----------------------------------------------------------------------
423      !!                  ***  ROUTINE dta_dyn  ***
424      !!
425      !! ** Purpose :  Prepares dynamics and physics fields from a NEMO run
426      !!               for an off-line simulation of passive tracers
427      !!
428      !! ** Method : calculates the position of data
429      !!             - computes slopes if needed
430      !!             - interpolates data if needed
431      !!----------------------------------------------------------------------
432      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
433      INTEGER, INTENT(in) ::   Kmm  ! ocean time level index
434      !
435      !!----------------------------------------------------------------------
436      !
437      IF( ln_timing )   CALL timing_start( 'dta_dyn_sed')
438      !
439      nsecdyn = nsec_year + nsec1jan000   ! number of seconds between Jan. 1st 00h of nit000 year and the middle of time step
440      !
441      IF( kt == nit000 ) THEN    ;    nprevrec = 0
442      ELSE                       ;    nprevrec = sf_dyn(jf_tem)%nrec_a(2)
443      ENDIF
444      CALL fld_read( kt, 1, sf_dyn )      !=  read data at kt time step   ==!
445      !
446      ts(:,:,:,jp_tem,Kmm) = sf_dyn(jf_tem)%fnow(:,:,:)  * tmask(:,:,:)    ! temperature
447      ts(:,:,:,jp_sal,Kmm) = sf_dyn(jf_sal)%fnow(:,:,:)  * tmask(:,:,:)    ! salinity
448      !
449      CALL eos    ( ts(:,:,:,:,Kmm), rhd, rhop, gdept_0(:,:,:) ) ! In any case, we need rhop
450
451      IF(ln_ctl) THEN                  ! print control
452         CALL prt_ctl(tab3d_1=ts(:,:,:,jp_tem,Kmm), clinfo1=' tn      - : ', mask1=tmask,  kdim=jpk   )
453         CALL prt_ctl(tab3d_1=ts(:,:,:,jp_sal,Kmm), clinfo1=' sn      - : ', mask1=tmask,  kdim=jpk   )
454      ENDIF
455      !
456      IF( ln_timing )   CALL timing_stop( 'dta_dyn_sed')
457      !
458   END SUBROUTINE dta_dyn_sed
459
460
461   SUBROUTINE dta_dyn_sed_init( Kmm )
462      !!----------------------------------------------------------------------
463      !!                  ***  ROUTINE dta_dyn_init  ***
464      !!
465      !! ** Purpose :   Initialisation of the dynamical data
466      !! ** Method  : - read the data namdta_dyn namelist
467      !!----------------------------------------------------------------------
468      INTEGER, INTENT( in ) :: Kmm                   ! ocean time level index
469      !
470      INTEGER  :: ierr, ierr0, ierr1, ierr2, ierr3   ! return error code
471      INTEGER  :: ifpr                               ! dummy loop indice
472      INTEGER  :: jfld                               ! dummy loop arguments
473      INTEGER  :: inum, idv, idimv                   ! local integer
474      INTEGER  :: ios                                ! Local integer output status for namelist read
475      !!
476      CHARACTER(len=100)            ::  cn_dir        !   Root directory for location of core files
477      TYPE(FLD_N), DIMENSION(2) ::  slf_d         ! array of namelist informations on the fields to read
478      TYPE(FLD_N) :: sn_tem , sn_sal   !   "                 "
479      !!
480      NAMELIST/namdta_dyn/cn_dir, ln_dynrnf, ln_dynrnf_depth,  fwbcorr, &
481         &                sn_tem, sn_sal
482      !!----------------------------------------------------------------------
483      !
484      REWIND( numnam_ref )              ! Namelist namdta_dyn in reference namelist : Offline: init. of dynamical data
485      READ  ( numnam_ref, namdta_dyn, IOSTAT = ios, ERR = 901)
486901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namdta_dyn in reference namelist', lwp )
487      REWIND( numnam_cfg )              ! Namelist namdta_dyn in configuration namelist : Offline: init. of dynamical data
488      READ  ( numnam_cfg, namdta_dyn, IOSTAT = ios, ERR = 902 )
489902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namdta_dyn in configuration namelist', lwp )
490      IF(lwm) WRITE ( numond, namdta_dyn )
491      !                                         ! store namelist information in an array
492      !                                         ! Control print
493      IF(lwp) THEN
494         WRITE(numout,*)
495         WRITE(numout,*) 'dta_dyn : offline dynamics '
496         WRITE(numout,*) '~~~~~~~ '
497         WRITE(numout,*) '   Namelist namdta_dyn'
498         WRITE(numout,*) '      runoffs option enabled (T) or not (F)            ln_dynrnf        = ', ln_dynrnf
499         WRITE(numout,*) '      runoffs is spread in vertical                    ln_dynrnf_depth  = ', ln_dynrnf_depth
500         WRITE(numout,*) '      annual global mean of empmr for ssh correction   fwbcorr          = ', fwbcorr
501         WRITE(numout,*)
502      ENDIF
503      !
504      jf_tem  = 1     ;   jf_sal  = 2    ;   jfld   = jf_sal
505      !
506      slf_d(jf_tem)  = sn_tem    ;   slf_d(jf_sal)  = sn_sal
507      !
508      ALLOCATE( sf_dyn(jfld), STAT=ierr )         ! set sf structure
509      IF( ierr > 0 )  THEN
510         CALL ctl_stop( 'dta_dyn: unable to allocate sf structure' )   ;   RETURN
511      ENDIF
512      !                                         ! fill sf with slf_i and control print
513      CALL fld_fill( sf_dyn, slf_d, cn_dir, 'dta_dyn_init', 'Data in file', 'namdta_dyn' )
514      !
515      ! Open file for each variable to get his number of dimension
516      DO ifpr = 1, jfld
517         CALL fld_clopn( sf_dyn(ifpr), nyear, nmonth, nday )
518         idv   = iom_varid( sf_dyn(ifpr)%num , slf_d(ifpr)%clvar )        ! id of the variable sdjf%clvar
519         idimv = iom_file ( sf_dyn(ifpr)%num )%ndims(idv)                 ! number of dimension for variable sdjf%clvar
520         IF( sf_dyn(ifpr)%num /= 0 )   CALL iom_close( sf_dyn(ifpr)%num ) ! close file if already open
521         ierr1=0
522         IF( idimv == 3 ) THEN    ! 2D variable
523                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,1)    , STAT=ierr0 )
524            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,1,2)  , STAT=ierr1 )
525         ELSE                     ! 3D variable
526                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,jpk)  , STAT=ierr0 )
527            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,jpk,2), STAT=ierr1 )
528         ENDIF
529         IF( ierr0 + ierr1 > 0 ) THEN
530            CALL ctl_stop( 'dta_dyn_init : unable to allocate sf_dyn array structure' )   ;   RETURN
531         ENDIF
532      END DO
533      !
534      CALL dta_dyn_sed( nit000, Kmm )
535      !
536   END SUBROUTINE dta_dyn_sed_init
537
538   SUBROUTINE dta_dyn_atf( kt, Kbb, Kmm, Kaa )
539     !!---------------------------------------------------------------------
540      !!                    ***  ROUTINE dta_dyn_swp  ***
541      !!
542      !! ** Purpose :   Asselin time filter of now SSH
543      !!---------------------------------------------------------------------
544      INTEGER, INTENT(in) :: kt             ! time step
545      INTEGER, INTENT(in) :: Kbb, Kmm, Kaa  ! ocean time level indices
546      !
547      !!---------------------------------------------------------------------
548
549      IF( kt == nit000 ) THEN
550         IF(lwp) WRITE(numout,*)
551         IF(lwp) WRITE(numout,*) 'dta_dyn_atf : Asselin time filter of sea surface height'
552         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ '
553      ENDIF
554
555      ssh(:,:,Kmm) = ssh(:,:,Kmm) + atfp * ( ssh(:,:,Kbb) - 2 * ssh(:,:,Kmm) + ssh(:,:,Kaa)) 
556
557      !! Do we also need to time filter e3t??
558      !
559   END SUBROUTINE dta_dyn_atf
560   
561   SUBROUTINE dta_dyn_sf_interp( kt, Kmm )
562      !!---------------------------------------------------------------------
563      !!                    ***  ROUTINE dta_dyn_sf_interp  ***
564      !!
565      !! ** Purpose :   Calculate scale factors at U/V/W points and depths
566      !!                given the after e3t field
567      !!---------------------------------------------------------------------
568      INTEGER, INTENT(in) :: kt   ! time step
569      INTEGER, INTENT(in) :: Kmm  ! ocean time level indices
570      !
571      INTEGER             :: ji, jj, jk
572      REAL(wp)            :: zcoef
573      !!---------------------------------------------------------------------
574
575      ! Horizontal scale factor interpolations
576      ! --------------------------------------
577      CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3u(:,:,:,Kmm), 'U' )
578      CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3v(:,:,:,Kmm), 'V' )
579
580      ! Vertical scale factor interpolations
581      ! ------------------------------------
582      CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3w (:,:,:,Kmm), 'W' )
583
584      ! t- and w- points depth
585      ! ----------------------
586      gdept(:,:,1,Kmm) = 0.5_wp * e3w(:,:,1,Kmm)
587      gdepw(:,:,1,Kmm) = 0.0_wp
588      !
589      DO jk = 2, jpk
590         DO jj = 1,jpj
591            DO ji = 1,jpi
592               zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
593               gdepw(ji,jj,jk,Kmm) = gdepw(ji,jj,jk-1,Kmm) + e3t(ji,jj,jk-1,Kmm)
594               gdept(ji,jj,jk,Kmm) =      zcoef  * ( gdepw(ji,jj,jk  ,Kmm) + 0.5 * e3w(ji,jj,jk,Kmm))  &
595                  &                + (1-zcoef) * ( gdept(ji,jj,jk-1,Kmm) + e3w(ji,jj,jk,Kmm))
596            END DO
597         END DO
598      END DO
599      !
600   END SUBROUTINE dta_dyn_sf_interp
601
602   SUBROUTINE dta_dyn_ssh( kt, phdivtr, psshb,  pemp, pssha, pe3ta )
603      !!----------------------------------------------------------------------
604      !!                ***  ROUTINE dta_dyn_wzv  ***
605      !!                   
606      !! ** Purpose :   compute the after ssh (ssh(:,:,Kaa)) and the now vertical velocity
607      !!
608      !! ** Method  : Using the incompressibility hypothesis,
609      !!        - the ssh increment is computed by integrating the horizontal divergence
610      !!          and multiply by the time step.
611      !!
612      !!        - compute the after scale factor : repartition of ssh INCREMENT proportionnaly
613      !!                                           to the level thickness ( z-star case )
614      !!
615      !!        - the vertical velocity is computed by integrating the horizontal divergence 
616      !!          from the bottom to the surface minus the scale factor evolution.
617      !!          The boundary conditions are w=0 at the bottom (no flux)
618      !!
619      !! ** action  :   ssh(:,:,Kaa) / e3t(:,:,:,Kaa) / ww
620      !!
621      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
622      !!----------------------------------------------------------------------
623      INTEGER,                                   INTENT(in )    :: kt        !  time-step
624      REAL(wp), DIMENSION(jpi,jpj,jpk)          , INTENT(in )   :: phdivtr   ! horizontal divergence transport
625      REAL(wp), DIMENSION(jpi,jpj)    , OPTIONAL, INTENT(in )   :: psshb     ! now ssh
626      REAL(wp), DIMENSION(jpi,jpj)    , OPTIONAL, INTENT(in )   :: pemp      ! evaporation minus precipitation
627      REAL(wp), DIMENSION(jpi,jpj)    , OPTIONAL, INTENT(inout) :: pssha     ! after ssh
628      REAL(wp), DIMENSION(jpi,jpj,jpk), OPTIONAL, INTENT(out)   :: pe3ta     ! after vertical scale factor
629      !
630      INTEGER                       :: jk
631      REAL(wp), DIMENSION(jpi,jpj)  :: zhdiv 
632      REAL(wp)  :: z2dt 
633      !!----------------------------------------------------------------------
634      !
635      z2dt = 2._wp * rdt
636      !
637      zhdiv(:,:) = 0._wp
638      DO jk = 1, jpkm1
639         zhdiv(:,:) = zhdiv(:,:) +  phdivtr(:,:,jk) * tmask(:,:,jk)
640      END DO
641      !                                                ! Sea surface  elevation time-stepping
642      pssha(:,:) = ( psshb(:,:) - z2dt * ( r1_rau0 * pemp(:,:)  + zhdiv(:,:) ) ) * ssmask(:,:)
643      !                                                 !
644      !                                                 ! After acale factors at t-points ( z_star coordinate )
645      DO jk = 1, jpkm1
646        pe3ta(:,:,jk) = e3t_0(:,:,jk) * ( 1._wp + pssha(:,:) * tmask(:,:,1) / ( ht_0(:,:) + 1.0 - tmask(:,:,1) ) )
647      END DO
648      !
649   END SUBROUTINE dta_dyn_ssh
650
651
652   SUBROUTINE dta_dyn_hrnf( Kmm )
653      !!----------------------------------------------------------------------
654      !!                  ***  ROUTINE sbc_rnf  ***
655      !!
656      !! ** Purpose :   update the horizontal divergence with the runoff inflow
657      !!
658      !! ** Method  :
659      !!                CAUTION : rnf is positive (inflow) decreasing the
660      !!                          divergence and expressed in m/s
661      !!
662      !! ** Action  :   phdivn   decreased by the runoff inflow
663      !!----------------------------------------------------------------------
664      !!
665      INTEGER, INTENT(in) :: Kmm  ! ocean time level index
666      !
667      INTEGER  ::   ji, jj, jk    ! dummy loop indices
668      !!----------------------------------------------------------------------
669      !
670      DO jj = 1, jpj                   ! update the depth over which runoffs are distributed
671         DO ji = 1, jpi
672            h_rnf(ji,jj) = 0._wp
673            DO jk = 1, nk_rnf(ji,jj)                           ! recalculates h_rnf to be the depth in metres
674                h_rnf(ji,jj) = h_rnf(ji,jj) + e3t(ji,jj,jk,Kmm)   ! to the bottom of the relevant grid box
675            END DO
676        END DO
677      END DO
678      !
679   END SUBROUTINE dta_dyn_hrnf
680
681
682
683   SUBROUTINE dta_dyn_slp( kt, Kbb, Kmm )
684      !!---------------------------------------------------------------------
685      !!                    ***  ROUTINE dta_dyn_slp  ***
686      !!
687      !! ** Purpose : Computation of slope
688      !!
689      !!---------------------------------------------------------------------
690      INTEGER,  INTENT(in) :: kt       ! time step
691      INTEGER,  INTENT(in) :: Kbb, Kmm ! ocean time level indices
692      !
693      INTEGER  ::   ji, jj     ! dummy loop indices
694      REAL(wp) ::   ztinta     ! ratio applied to after  records when doing time interpolation
695      REAL(wp) ::   ztintb     ! ratio applied to before records when doing time interpolation
696      INTEGER  ::   iswap 
697      REAL(wp), DIMENSION(jpi,jpj,jpk)      ::   zuslp, zvslp, zwslpi, zwslpj
698      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts) ::   zts
699      !!---------------------------------------------------------------------
700      !
701      IF( sf_dyn(jf_tem)%ln_tint ) THEN    ! Computes slopes (here avt is used as workspace)                       
702         IF( kt == nit000 ) THEN
703            IF(lwp) WRITE(numout,*) ' Compute new slopes at kt = ', kt
704            zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,1) * tmask(:,:,:)   ! temperature
705            zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,1) * tmask(:,:,:)   ! salinity
706            avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,1) * tmask(:,:,:)   ! vertical diffusive coef.
707            CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
708            uslpdta (:,:,:,1) = zuslp (:,:,:) 
709            vslpdta (:,:,:,1) = zvslp (:,:,:) 
710            wslpidta(:,:,:,1) = zwslpi(:,:,:) 
711            wslpjdta(:,:,:,1) = zwslpj(:,:,:) 
712            !
713            zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,2) * tmask(:,:,:)   ! temperature
714            zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,2) * tmask(:,:,:)   ! salinity
715            avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,2) * tmask(:,:,:)   ! vertical diffusive coef.
716            CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
717            uslpdta (:,:,:,2) = zuslp (:,:,:) 
718            vslpdta (:,:,:,2) = zvslp (:,:,:) 
719            wslpidta(:,:,:,2) = zwslpi(:,:,:) 
720            wslpjdta(:,:,:,2) = zwslpj(:,:,:) 
721         ELSE
722           !
723           iswap = 0
724           IF( sf_dyn(jf_tem)%nrec_a(2) - nprevrec /= 0 )  iswap = 1
725           IF( nsecdyn > sf_dyn(jf_tem)%nrec_b(2) .AND. iswap == 1 )  THEN    ! read/update the after data
726              IF(lwp) WRITE(numout,*) ' Compute new slopes at kt = ', kt
727              uslpdta (:,:,:,1) =  uslpdta (:,:,:,2)         ! swap the data
728              vslpdta (:,:,:,1) =  vslpdta (:,:,:,2) 
729              wslpidta(:,:,:,1) =  wslpidta(:,:,:,2) 
730              wslpjdta(:,:,:,1) =  wslpjdta(:,:,:,2) 
731              !
732              zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,2) * tmask(:,:,:)   ! temperature
733              zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,2) * tmask(:,:,:)   ! salinity
734              avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,2) * tmask(:,:,:)   ! vertical diffusive coef.
735              CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
736              !
737              uslpdta (:,:,:,2) = zuslp (:,:,:) 
738              vslpdta (:,:,:,2) = zvslp (:,:,:) 
739              wslpidta(:,:,:,2) = zwslpi(:,:,:) 
740              wslpjdta(:,:,:,2) = zwslpj(:,:,:) 
741            ENDIF
742         ENDIF
743      ENDIF
744      !
745      IF( sf_dyn(jf_tem)%ln_tint )  THEN
746         ztinta =  REAL( nsecdyn - sf_dyn(jf_tem)%nrec_b(2), wp )  &
747            &    / REAL( sf_dyn(jf_tem)%nrec_a(2) - sf_dyn(jf_tem)%nrec_b(2), wp )
748         ztintb =  1. - ztinta
749         IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN    ! Computes slopes (here avt is used as workspace)
750            uslp (:,:,:) = ztintb * uslpdta (:,:,:,1)  + ztinta * uslpdta (:,:,:,2) 
751            vslp (:,:,:) = ztintb * vslpdta (:,:,:,1)  + ztinta * vslpdta (:,:,:,2) 
752            wslpi(:,:,:) = ztintb * wslpidta(:,:,:,1)  + ztinta * wslpidta(:,:,:,2) 
753            wslpj(:,:,:) = ztintb * wslpjdta(:,:,:,1)  + ztinta * wslpjdta(:,:,:,2) 
754         ENDIF
755      ELSE
756         zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fnow(:,:,:) * tmask(:,:,:)   ! temperature
757         zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fnow(:,:,:) * tmask(:,:,:)   ! salinity
758         avt(:,:,:)        = sf_dyn(jf_avt)%fnow(:,:,:) * tmask(:,:,:)   ! vertical diffusive coef.
759         CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
760         !
761         IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN    ! Computes slopes (here avt is used as workspace)
762            uslp (:,:,:) = zuslp (:,:,:)
763            vslp (:,:,:) = zvslp (:,:,:)
764            wslpi(:,:,:) = zwslpi(:,:,:)
765            wslpj(:,:,:) = zwslpj(:,:,:)
766         ENDIF
767      ENDIF
768      !
769   END SUBROUTINE dta_dyn_slp
770
771
772   SUBROUTINE compute_slopes( kt, pts, puslp, pvslp, pwslpi, pwslpj, Kbb, Kmm )
773      !!---------------------------------------------------------------------
774      !!                    ***  ROUTINE dta_dyn_slp  ***
775      !!
776      !! ** Purpose :   Computation of slope
777      !!---------------------------------------------------------------------
778      INTEGER ,                              INTENT(in ) :: kt       ! time step
779      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in ) :: pts      ! temperature/salinity
780      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: puslp    ! zonal isopycnal slopes
781      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pvslp    ! meridional isopycnal slopes
782      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pwslpi   ! zonal diapycnal slopes
783      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pwslpj   ! meridional diapycnal slopes
784      INTEGER ,                              INTENT(in ) :: Kbb, Kmm ! ocean time level indices
785      !!---------------------------------------------------------------------
786      !
787      IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN    ! Computes slopes (here avt is used as workspace)
788         CALL eos    ( pts, rhd, rhop, gdept_0(:,:,:) )
789         CALL eos_rab( pts, rab_n, Kmm )       ! now local thermal/haline expension ratio at T-points
790         CALL bn2    ( pts, rab_n, rn2, Kmm  ) ! now    Brunt-Vaisala
791
792      ! Partial steps: before Horizontal DErivative
793      IF( ln_zps  .AND. .NOT. ln_isfcav)                            &
794         &            CALL zps_hde    ( kt, Kmm, jpts, pts, gtsu, gtsv,  &  ! Partial steps: before horizontal gradient
795         &                                        rhd, gru , grv    )  ! of t, s, rd at the last ocean level
796      IF( ln_zps .AND.        ln_isfcav)                            &
797         &            CALL zps_hde_isf( kt, Kmm, jpts, pts, gtsu, gtsv, gtui, gtvi, &  ! Partial steps for top cell (ISF)
798         &                                        rhd, gru , grv , grui, grvi )  ! of t, s, rd at the first ocean level
799
800         rn2b(:,:,:) = rn2(:,:,:)                ! needed for zdfmxl
801         CALL zdf_mxl( kt, Kmm )                 ! mixed layer depth
802         CALL ldf_slp( kt, rhd, rn2, Kbb, Kmm )  ! slopes
803         puslp (:,:,:) = uslp (:,:,:)
804         pvslp (:,:,:) = vslp (:,:,:)
805         pwslpi(:,:,:) = wslpi(:,:,:)
806         pwslpj(:,:,:) = wslpj(:,:,:)
807     ELSE
808         puslp (:,:,:) = 0.            ! to avoid warning when compiling
809         pvslp (:,:,:) = 0.
810         pwslpi(:,:,:) = 0.
811         pwslpj(:,:,:) = 0.
812     ENDIF
813      !
814   END SUBROUTINE compute_slopes
815
816   !!======================================================================
817END MODULE dtadyn
Note: See TracBrowser for help on using the repository browser.