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

source: NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/OFF/dtadyn.F90 @ 12246

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

rev12232_dev_r12072_MERGE_OPTION2_2019: add modifications from dev_r12114_ticket_2263, results unchanged except SPITZ12 as explained in #2263

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