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

source: NEMO/trunk/src/OFF/dtadyn.F90 @ 13295

Last change on this file since 13295 was 13295, checked in by acc, 4 years ago

Replace do-loop macros in the trunk with alternative forms with greater flexibility for extra halo applications. This alters a lot of routines but does not change any behaviour or results. do_loop_substitute.h90 is greatly simplified by this change. SETTE results are identical to those with the previous revision

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