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_r11943_MERGE_2019/src/OFF – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OFF/dtadyn.F90 @ 12250

Last change on this file since 12250 was 12250, checked in by smasson, 4 years ago

rev12240_dev_r11943_MERGE_2019: same as [12246], add modifications from dev_r12114_ticket_2263, results unchanged except SPITZ12 as explained in #2263

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