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

source: NEMO/branches/2019/dev_r11351_fldread_with_XIOS/src/OFF/dtadyn.F90 @ 13463

Last change on this file since 13463 was 13463, checked in by andmirek, 4 years ago

Ticket #2195:update to trunk 13461

  • Property svn:keywords set to Id
File size: 42.8 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      !
299      IF( ln_trabbl ) THEN
300         sf_dyn(jf_ubl)%cltype = 'U'   ;   sf_dyn(jf_ubl)%zsgn =  1._wp 
301         sf_dyn(jf_vbl)%cltype = 'V'   ;   sf_dyn(jf_vbl)%zsgn =  1._wp 
302      END IF
303      !
304      ! Open file for each variable to get his number of dimension
305      DO ifpr = 1, jfld
306         CALL fld_def( sf_dyn(ifpr) )
307         CALL iom_open( sf_dyn(ifpr)%clname, sf_dyn(ifpr)%num )
308         idv   = iom_varid( sf_dyn(ifpr)%num , slf_d(ifpr)%clvar )        ! id of the variable sdjf%clvar
309         idimv = iom_file ( sf_dyn(ifpr)%num )%ndims(idv)                 ! number of dimension for variable sdjf%clvar
310         CALL iom_close( sf_dyn(ifpr)%num )                               ! close file if already open
311         ierr1=0
312         IF( idimv == 3 ) THEN    ! 2D variable
313                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,1)    , STAT=ierr0 )
314            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,1,2)  , STAT=ierr1 )
315         ELSE                     ! 3D variable
316                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,jpk)  , STAT=ierr0 )
317            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,jpk,2), STAT=ierr1 )
318         ENDIF
319         IF( ierr0 + ierr1 > 0 ) THEN
320            CALL ctl_stop( 'dta_dyn_init : unable to allocate sf_dyn array structure' )   ;   RETURN
321         ENDIF
322      END DO
323      !
324      IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN                  ! slopes
325         IF( sf_dyn(jf_tem)%ln_tint ) THEN      ! time interpolation
326            ALLOCATE( uslpdta (jpi,jpj,jpk,2), vslpdta (jpi,jpj,jpk,2),    &
327            &         wslpidta(jpi,jpj,jpk,2), wslpjdta(jpi,jpj,jpk,2), STAT=ierr2 )
328            !
329            IF( ierr2 > 0 )  THEN
330               CALL ctl_stop( 'dta_dyn_init : unable to allocate slope arrays' )   ;   RETURN
331            ENDIF
332         ENDIF
333      ENDIF
334      !
335      IF( .NOT.ln_linssh ) THEN
336        IF( .NOT. sf_dyn(jf_uwd)%ln_clim .AND. ln_rsttr .AND.    &                     ! Restart: read in restart file
337           iom_varid( numrtr, 'sshn', ldstop = .FALSE. ) > 0 ) THEN
338           IF(lwp) WRITE(numout,*) ' ssh(:,:,Kmm) forcing fields read in the restart file for initialisation'
339           CALL iom_get( numrtr, jpdom_auto, 'sshn', ssh(:,:,Kmm)   )
340           CALL iom_get( numrtr, jpdom_auto, 'sshb', ssh(:,:,Kbb)   )
341        ELSE
342           IF(lwp) WRITE(numout,*) ' ssh(:,:,Kmm) forcing fields read in the restart file for initialisation'
343           CALL iom_open( 'restart', inum )
344           CALL iom_get( inum, jpdom_auto, 'sshn', ssh(:,:,Kmm)   )
345           CALL iom_get( inum, jpdom_auto, 'sshb', ssh(:,:,Kbb)   )
346           CALL iom_close( inum )                                        ! close file
347        ENDIF
348        !
349#if defined key_qco
350        CALL dom_qco_r3c( ssh(:,:,Kbb), r3t(:,:,Kbb), r3u(:,:,Kbb), r3v(:,:,Kbb) )
351        CALL dom_qco_r3c( ssh(:,:,Kmm), r3t(:,:,Kmm), r3u(:,:,Kmm), r3v(:,:,Kmm) )
352#else
353        DO jk = 1, jpkm1
354           e3t(:,:,jk,Kmm) = e3t_0(:,:,jk) * ( 1._wp + ssh(:,:,Kmm) * r1_ht_0(:,:) * tmask(:,:,jk) )
355        ENDDO
356        e3t(:,:,jpk,Kaa) = e3t_0(:,:,jpk)
357
358        ! Horizontal scale factor interpolations
359        ! --------------------------------------
360        CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3u(:,:,:,Kmm), 'U' )
361        CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3v(:,:,:,Kmm), 'V' )
362
363        ! Vertical scale factor interpolations
364        ! ------------------------------------
365        CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3w(:,:,:,Kmm), 'W' )
366!!gm this should be computed from ssh(Kbb) 
367        e3t(:,:,:,Kbb)  = e3t(:,:,:,Kmm)
368        e3u(:,:,:,Kbb)  = e3u(:,:,:,Kmm)
369        e3v(:,:,:,Kbb)  = e3v(:,:,:,Kmm)
370
371        ! t- and w- points depth
372        ! ----------------------
373        gdept(:,:,1,Kmm) = 0.5_wp * e3w(:,:,1,Kmm)
374        gdepw(:,:,1,Kmm) = 0.0_wp
375
376        DO_3D( 1, 1, 1, 1, 2, jpk )
377          !    zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))   ! 0 everywhere
378          !    tmask = wmask, ie everywhere expect at jk = mikt
379                                                             ! 1 for jk =
380                                                             ! mikt
381           zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
382           gdepw(ji,jj,jk,Kmm) = gdepw(ji,jj,jk-1,Kmm) + e3t(ji,jj,jk-1,Kmm)
383           gdept(ji,jj,jk,Kmm) =      zcoef  * ( gdepw(ji,jj,jk  ,Kmm) + 0.5 * e3w(ji,jj,jk,Kmm))  &
384               &                + (1-zcoef) * ( gdept(ji,jj,jk-1,Kmm) + e3w(ji,jj,jk,Kmm))
385        END_3D
386
387        gdept(:,:,:,Kbb) = gdept(:,:,:,Kmm)
388        gdepw(:,:,:,Kbb) = gdepw(:,:,:,Kmm)
389        !
390      ENDIF
391#endif
392      !
393      IF( ln_dynrnf .AND. ln_dynrnf_depth ) THEN       ! read depht over which runoffs are distributed
394         IF(lwp) WRITE(numout,*) 
395         IF(lwp) WRITE(numout,*) ' read in the file depht over which runoffs are distributed'
396         CALL iom_open ( "runoffs", inum )                           ! open file
397         CALL iom_get  ( inum, jpdom_global, 'rodepth', h_rnf )   ! read the river mouth array
398         CALL iom_close( inum )                                        ! close file
399         !
400         nk_rnf(:,:) = 0                               ! set the number of level over which river runoffs are applied
401         DO_2D( 1, 1, 1, 1 )
402            IF( h_rnf(ji,jj) > 0._wp ) THEN
403               jk = 2
404               DO WHILE ( jk /= mbkt(ji,jj) .AND. gdept_0(ji,jj,jk) < h_rnf(ji,jj) ) ;  jk = jk + 1
405               END DO
406               nk_rnf(ji,jj) = jk
407            ELSEIF( h_rnf(ji,jj) == -1._wp   ) THEN   ;  nk_rnf(ji,jj) = 1
408            ELSEIF( h_rnf(ji,jj) == -999._wp ) THEN   ;  nk_rnf(ji,jj) = mbkt(ji,jj)
409            ELSE
410               CALL ctl_stop( 'sbc_rnf_init: runoff depth not positive, and not -999 or -1, rnf value in file fort.999'  )
411               WRITE(999,*) 'ji, jj, h_rnf(ji,jj) :', ji, jj, h_rnf(ji,jj)
412            ENDIF
413         END_2D
414!!st pourquoi on n'utilise pas le gde3w ici plutôt que de faire une boucle ?
415         DO_2D( 1, 1, 1, 1 )
416            h_rnf(ji,jj) = 0._wp
417            DO jk = 1, nk_rnf(ji,jj)
418               h_rnf(ji,jj) = h_rnf(ji,jj) + e3t(ji,jj,jk,Kmm)
419            END DO
420         END_2D
421      ELSE                                       ! runoffs applied at the surface
422         nk_rnf(:,:) = 1
423         h_rnf (:,:) = e3t(:,:,1,Kmm)
424      ENDIF
425      nkrnf_max = MAXVAL( nk_rnf(:,:) )
426      hrnf_max = MAXVAL( h_rnf(:,:) )
427      IF( lk_mpp )  THEN
428         CALL mpp_max( 'dtadyn', nkrnf_max )                 ! max over the  global domain
429         CALL mpp_max( 'dtadyn', hrnf_max )                 ! max over the  global domain
430      ENDIF
431      IF(lwp) WRITE(numout,*) ' '
432      IF(lwp) WRITE(numout,*) ' max depht of runoff : ', hrnf_max,'    max level  : ', nkrnf_max
433      IF(lwp) WRITE(numout,*) ' '
434      !
435      CALL dta_dyn( nit000, Kbb, Kmm, Kaa )
436      !
437   END SUBROUTINE dta_dyn_init
438
439   
440   SUBROUTINE dta_dyn_sed( kt, Kmm )
441      !!----------------------------------------------------------------------
442      !!                  ***  ROUTINE dta_dyn  ***
443      !!
444      !! ** Purpose :  Prepares dynamics and physics fields from a NEMO run
445      !!               for an off-line simulation of passive tracers
446      !!
447      !! ** Method : calculates the position of data
448      !!             - computes slopes if needed
449      !!             - interpolates data if needed
450      !!----------------------------------------------------------------------
451      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
452      INTEGER, INTENT(in) ::   Kmm  ! ocean time level index
453      !
454      !!----------------------------------------------------------------------
455      !
456      IF( ln_timing )   CALL timing_start( 'dta_dyn_sed')
457      !
458      nsecdyn = nsec_year + nsec1jan000   ! number of seconds between Jan. 1st 00h of nit000 year and the middle of time step
459      !
460      IF( kt == nit000 ) THEN    ;    nprevrec = 0
461      ELSE                       ;    nprevrec = sf_dyn(jf_tem)%nrec(2,sf_dyn(jf_tem)%naa)
462      ENDIF
463      CALL fld_read( kt, 1, sf_dyn )      !=  read data at kt time step   ==!
464      !
465      ts(:,:,:,jp_tem,Kmm) = sf_dyn(jf_tem)%fnow(:,:,:)  * tmask(:,:,:)    ! temperature
466      ts(:,:,:,jp_sal,Kmm) = sf_dyn(jf_sal)%fnow(:,:,:)  * tmask(:,:,:)    ! salinity
467      !
468      CALL eos    ( ts(:,:,:,:,Kmm), rhd, rhop, gdept_0(:,:,:) ) ! In any case, we need rhop
469
470      IF(sn_cfctl%l_prtctl) THEN                     ! print control
471         CALL prt_ctl(tab3d_1=ts(:,:,:,jp_tem,Kmm), clinfo1=' tn      - : ', mask1=tmask,  kdim=jpk   )
472         CALL prt_ctl(tab3d_1=ts(:,:,:,jp_sal,Kmm), clinfo1=' sn      - : ', mask1=tmask,  kdim=jpk   )
473      ENDIF
474      !
475      IF( ln_timing )   CALL timing_stop( 'dta_dyn_sed')
476      !
477   END SUBROUTINE dta_dyn_sed
478
479
480   SUBROUTINE dta_dyn_sed_init( Kmm )
481      !!----------------------------------------------------------------------
482      !!                  ***  ROUTINE dta_dyn_init  ***
483      !!
484      !! ** Purpose :   Initialisation of the dynamical data
485      !! ** Method  : - read the data namdta_dyn namelist
486      !!----------------------------------------------------------------------
487      INTEGER, INTENT( in ) :: Kmm                   ! ocean time level index
488      !
489      INTEGER  :: ierr, ierr0, ierr1, ierr2, ierr3   ! return error code
490      INTEGER  :: ifpr                               ! dummy loop indice
491      INTEGER  :: jfld                               ! dummy loop arguments
492      INTEGER  :: inum, idv, idimv                   ! local integer
493      INTEGER  :: ios                                ! Local integer output status for namelist read
494      !!
495      CHARACTER(len=100)            ::  cn_dir        !   Root directory for location of core files
496      TYPE(FLD_N), DIMENSION(2) ::  slf_d         ! array of namelist informations on the fields to read
497      TYPE(FLD_N) :: sn_tem , sn_sal   !   "                 "
498      !!
499      NAMELIST/namdta_dyn/cn_dir, ln_dynrnf, ln_dynrnf_depth,  fwbcorr, &
500         &                sn_tem, sn_sal
501      !!----------------------------------------------------------------------
502      !
503      READ  ( numnam_ref, namdta_dyn, IOSTAT = ios, ERR = 901)
504901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namdta_dyn in reference namelist' )
505      READ  ( numnam_cfg, namdta_dyn, IOSTAT = ios, ERR = 902 )
506902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namdta_dyn in configuration namelist' )
507      IF(lwm) WRITE ( numond, namdta_dyn )
508      !                                         ! store namelist information in an array
509      !                                         ! Control print
510      IF(lwp) THEN
511         WRITE(numout,*)
512         WRITE(numout,*) 'dta_dyn : offline dynamics '
513         WRITE(numout,*) '~~~~~~~ '
514         WRITE(numout,*) '   Namelist namdta_dyn'
515         WRITE(numout,*) '      runoffs option enabled (T) or not (F)            ln_dynrnf        = ', ln_dynrnf
516         WRITE(numout,*) '      runoffs is spread in vertical                    ln_dynrnf_depth  = ', ln_dynrnf_depth
517         WRITE(numout,*) '      annual global mean of empmr for ssh correction   fwbcorr          = ', fwbcorr
518         WRITE(numout,*)
519      ENDIF
520      !
521      jf_tem  = 1     ;   jf_sal  = 2    ;   jfld   = jf_sal
522      !
523      slf_d(jf_tem)  = sn_tem    ;   slf_d(jf_sal)  = sn_sal
524      !
525      ALLOCATE( sf_dyn(jfld), STAT=ierr )         ! set sf structure
526      IF( ierr > 0 )  THEN
527         CALL ctl_stop( 'dta_dyn: unable to allocate sf structure' )   ;   RETURN
528      ENDIF
529      !                                         ! fill sf with slf_i and control print
530      CALL fld_fill( sf_dyn, slf_d, cn_dir, 'dta_dyn_init', 'Data in file', 'namdta_dyn' )
531      !
532      ! Open file for each variable to get his number of dimension
533      DO ifpr = 1, jfld
534         CALL fld_def( sf_dyn(ifpr) )
535         CALL iom_open( sf_dyn(ifpr)%clname, sf_dyn(ifpr)%num )
536         idv   = iom_varid( sf_dyn(ifpr)%num , slf_d(ifpr)%clvar )        ! id of the variable sdjf%clvar
537         idimv = iom_file ( sf_dyn(ifpr)%num )%ndims(idv)                 ! number of dimension for variable sdjf%clvar
538         CALL iom_close( sf_dyn(ifpr)%num )                               ! close file if already open
539         ierr1=0
540         IF( idimv == 3 ) THEN    ! 2D variable
541                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,1)    , STAT=ierr0 )
542            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,1,2)  , STAT=ierr1 )
543         ELSE                     ! 3D variable
544                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,jpk)  , STAT=ierr0 )
545            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,jpk,2), STAT=ierr1 )
546         ENDIF
547         IF( ierr0 + ierr1 > 0 ) THEN
548            CALL ctl_stop( 'dta_dyn_init : unable to allocate sf_dyn array structure' )   ;   RETURN
549         ENDIF
550      END DO
551      !
552      CALL dta_dyn_sed( nit000, Kmm )
553      !
554   END SUBROUTINE dta_dyn_sed_init
555
556   
557   SUBROUTINE dta_dyn_atf( kt, Kbb, Kmm, Kaa )
558     !!---------------------------------------------------------------------
559      !!                    ***  ROUTINE dta_dyn_swp  ***
560      !!
561      !! ** Purpose :   Asselin time filter of now SSH
562      !!---------------------------------------------------------------------
563      INTEGER, INTENT(in) :: kt             ! time step
564      INTEGER, INTENT(in) :: Kbb, Kmm, Kaa  ! ocean time level indices
565      !
566      !!---------------------------------------------------------------------
567
568      IF( kt == nit000 ) THEN
569         IF(lwp) WRITE(numout,*)
570         IF(lwp) WRITE(numout,*) 'dta_dyn_atf : Asselin time filter of sea surface height'
571         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ '
572      ENDIF
573
574      ssh(:,:,Kmm) = ssh(:,:,Kmm) + rn_atfp * ( ssh(:,:,Kbb) - 2 * ssh(:,:,Kmm) + ssh(:,:,Kaa)) 
575
576      !! Do we also need to time filter e3t??
577      !
578   END SUBROUTINE dta_dyn_atf
579   
580   
581#if ! defined key_qco   
582   SUBROUTINE dta_dyn_sf_interp( kt, Kmm )
583      !!---------------------------------------------------------------------
584      !!                    ***  ROUTINE dta_dyn_sf_interp  ***
585      !!
586      !! ** Purpose :   Calculate scale factors at U/V/W points and depths
587      !!                given the after e3t field
588      !!---------------------------------------------------------------------
589      INTEGER, INTENT(in) :: kt   ! time step
590      INTEGER, INTENT(in) :: Kmm  ! ocean time level indices
591      !
592      INTEGER             :: ji, jj, jk
593      REAL(wp)            :: zcoef
594      !!---------------------------------------------------------------------
595
596      ! Horizontal scale factor interpolations
597      ! --------------------------------------
598      CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3u(:,:,:,Kmm), 'U' )
599      CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3v(:,:,:,Kmm), 'V' )
600
601      ! Vertical scale factor interpolations
602      ! ------------------------------------
603      CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3w (:,:,:,Kmm), 'W' )
604
605      ! t- and w- points depth
606      ! ----------------------
607      gdept(:,:,1,Kmm) = 0.5_wp * e3w(:,:,1,Kmm)
608      gdepw(:,:,1,Kmm) = 0.0_wp
609      !
610      DO_3D( 1, 1, 1, 1, 2, jpk )
611         zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
612         gdepw(ji,jj,jk,Kmm) = gdepw(ji,jj,jk-1,Kmm) + e3t(ji,jj,jk-1,Kmm)
613         gdept(ji,jj,jk,Kmm) =      zcoef  * ( gdepw(ji,jj,jk  ,Kmm) + 0.5 * e3w(ji,jj,jk,Kmm))  &
614            &                + (1-zcoef) * ( gdept(ji,jj,jk-1,Kmm) + e3w(ji,jj,jk,Kmm))
615      END_3D
616      !
617   END SUBROUTINE dta_dyn_sf_interp
618#endif
619
620
621   SUBROUTINE dta_dyn_ssh( kt, phdivtr, psshb,  pemp, pssha, pe3ta )
622      !!----------------------------------------------------------------------
623      !!                ***  ROUTINE dta_dyn_wzv  ***
624      !!                   
625      !! ** Purpose :   compute the after ssh (ssh(:,:,Kaa)) and the now vertical velocity
626      !!
627      !! ** Method  : Using the incompressibility hypothesis,
628      !!        - the ssh increment is computed by integrating the horizontal divergence
629      !!          and multiply by the time step.
630      !!
631      !!        - compute the after scale factor : repartition of ssh INCREMENT proportionnaly
632      !!                                           to the level thickness ( z-star case )
633      !!
634      !!        - the vertical velocity is computed by integrating the horizontal divergence 
635      !!          from the bottom to the surface minus the scale factor evolution.
636      !!          The boundary conditions are w=0 at the bottom (no flux)
637      !!
638      !! ** action  :   ssh(:,:,Kaa) / e3t(:,:,k,Kaa) / ww
639      !!
640      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
641      !!----------------------------------------------------------------------
642      INTEGER,                                   INTENT(in )    :: kt        !  time-step
643      REAL(wp), DIMENSION(jpi,jpj,jpk)          , INTENT(in )   :: phdivtr   ! horizontal divergence transport
644      REAL(wp), DIMENSION(jpi,jpj)    , OPTIONAL, INTENT(in )   :: psshb     ! now ssh
645      REAL(wp), DIMENSION(jpi,jpj)    , OPTIONAL, INTENT(in )   :: pemp      ! evaporation minus precipitation
646      REAL(wp), DIMENSION(jpi,jpj)    , OPTIONAL, INTENT(inout) :: pssha     ! after ssh
647      REAL(wp), DIMENSION(jpi,jpj,jpk), OPTIONAL, INTENT(out)   :: pe3ta     ! after vertical scale factor
648      !
649      INTEGER                       :: jk
650      REAL(wp), DIMENSION(jpi,jpj)  :: zhdiv 
651      REAL(wp)  :: z2dt 
652      !!----------------------------------------------------------------------
653      !
654      z2dt = 2._wp * rn_Dt
655      !
656      zhdiv(:,:) = 0._wp
657      DO jk = 1, jpkm1
658         zhdiv(:,:) = zhdiv(:,:) +  phdivtr(:,:,jk) * tmask(:,:,jk)
659      END DO
660      !                                                ! Sea surface  elevation time-stepping
661      pssha(:,:) = ( psshb(:,:) - z2dt * ( r1_rho0 * pemp(:,:)  + zhdiv(:,:) ) ) * ssmask(:,:)
662      !
663      IF( PRESENT( pe3ta ) ) THEN                      ! After acale factors at t-points ( z_star coordinate )
664      DO jk = 1, jpkm1
665            pe3ta(:,:,jk) = e3t_0(:,:,jk) * ( 1._wp + pssha(:,:) * r1_ht_0(:,:) * tmask(:,:,jk) )
666      END DO
667      ENDIF
668      !
669   END SUBROUTINE dta_dyn_ssh
670
671
672   SUBROUTINE dta_dyn_hrnf( Kmm )
673      !!----------------------------------------------------------------------
674      !!                  ***  ROUTINE sbc_rnf  ***
675      !!
676      !! ** Purpose :   update the horizontal divergence with the runoff inflow
677      !!
678      !! ** Method  :
679      !!                CAUTION : rnf is positive (inflow) decreasing the
680      !!                          divergence and expressed in m/s
681      !!
682      !! ** Action  :   phdivn   decreased by the runoff inflow
683      !!----------------------------------------------------------------------
684      !!
685      INTEGER, INTENT(in) :: Kmm  ! ocean time level index
686      !
687      INTEGER  ::   ji, jj, jk    ! dummy loop indices
688      !!----------------------------------------------------------------------
689      !
690!!st code dupliqué même remarque que plus haut pourquoi ne pas utiliser gdepw ?
691      DO_2D( 1, 1, 1, 1 )
692         h_rnf(ji,jj) = 0._wp
693         DO jk = 1, nk_rnf(ji,jj)                           ! recalculates h_rnf to be the depth in metres
694             h_rnf(ji,jj) = h_rnf(ji,jj) + e3t(ji,jj,jk,Kmm)   ! to the bottom of the relevant grid box
695         END DO
696      END_2D
697      !
698   END SUBROUTINE dta_dyn_hrnf
699
700
701
702   SUBROUTINE dta_dyn_slp( kt, Kbb, Kmm )
703      !!---------------------------------------------------------------------
704      !!                    ***  ROUTINE dta_dyn_slp  ***
705      !!
706      !! ** Purpose : Computation of slope
707      !!
708      !!---------------------------------------------------------------------
709      INTEGER,  INTENT(in) :: kt       ! time step
710      INTEGER,  INTENT(in) :: Kbb, Kmm ! ocean time level indices
711      !
712      INTEGER  ::   ji, jj     ! dummy loop indices
713      REAL(wp) ::   ztinta     ! ratio applied to after  records when doing time interpolation
714      REAL(wp) ::   ztintb     ! ratio applied to before records when doing time interpolation
715      INTEGER  ::   iswap 
716      REAL(wp), DIMENSION(jpi,jpj,jpk)      ::   zuslp, zvslp, zwslpi, zwslpj
717      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts) ::   zts
718      !!---------------------------------------------------------------------
719      !
720      IF( sf_dyn(jf_tem)%ln_tint ) THEN    ! Computes slopes (here avt is used as workspace)
721         !
722         IF( kt == nit000 ) THEN
723            IF(lwp) WRITE(numout,*) ' Compute new slopes at kt = ', kt
724            zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,sf_dyn(jf_tem)%nbb) * tmask(:,:,:)   ! temperature
725            zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,sf_dyn(jf_sal)%nbb) * tmask(:,:,:)   ! salinity
726            avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,sf_dyn(jf_avt)%nbb) * tmask(:,:,:)   ! vertical diffusive coef.
727            CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
728            uslpdta (:,:,:,1) = zuslp (:,:,:) 
729            vslpdta (:,:,:,1) = zvslp (:,:,:) 
730            wslpidta(:,:,:,1) = zwslpi(:,:,:) 
731            wslpjdta(:,:,:,1) = zwslpj(:,:,:) 
732            !
733            zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,sf_dyn(jf_tem)%naa) * tmask(:,:,:)   ! temperature
734            zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,sf_dyn(jf_sal)%naa) * tmask(:,:,:)   ! salinity
735            avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,sf_dyn(jf_avt)%naa) * tmask(:,:,:)   ! vertical diffusive coef.
736            CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
737            uslpdta (:,:,:,2) = zuslp (:,:,:) 
738            vslpdta (:,:,:,2) = zvslp (:,:,:) 
739            wslpidta(:,:,:,2) = zwslpi(:,:,:) 
740            wslpjdta(:,:,:,2) = zwslpj(:,:,:) 
741         ELSE
742           !
743           iswap = 0
744           IF( sf_dyn(jf_tem)%nrec(2,sf_dyn(jf_tem)%naa) - nprevrec /= 0 )  iswap = 1
745           IF( nsecdyn > sf_dyn(jf_tem)%nrec(2,sf_dyn(jf_tem)%nbb) .AND. iswap == 1 )  THEN    ! read/update the after data
746              IF(lwp) WRITE(numout,*) ' Compute new slopes at kt = ', kt
747              uslpdta (:,:,:,1) =  uslpdta (:,:,:,2)         ! swap the data
748              vslpdta (:,:,:,1) =  vslpdta (:,:,:,2) 
749              wslpidta(:,:,:,1) =  wslpidta(:,:,:,2) 
750              wslpjdta(:,:,:,1) =  wslpjdta(:,:,:,2) 
751              !
752              zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,sf_dyn(jf_tem)%naa) * tmask(:,:,:)   ! temperature
753              zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,sf_dyn(jf_sal)%naa) * tmask(:,:,:)   ! salinity
754              avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,sf_dyn(jf_avt)%naa) * tmask(:,:,:)   ! vertical diffusive coef.
755              CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
756              !
757              uslpdta (:,:,:,2) = zuslp (:,:,:) 
758              vslpdta (:,:,:,2) = zvslp (:,:,:) 
759              wslpidta(:,:,:,2) = zwslpi(:,:,:) 
760              wslpjdta(:,:,:,2) = zwslpj(:,:,:) 
761            ENDIF
762         ENDIF
763      ENDIF
764      !
765      IF( sf_dyn(jf_tem)%ln_tint )  THEN
766         ztinta =  REAL( nsecdyn - sf_dyn(jf_tem)%nrec(2,sf_dyn(jf_tem)%nbb), wp )  &
767            &    / REAL( sf_dyn(jf_tem)%nrec(2,sf_dyn(jf_tem)%naa) - sf_dyn(jf_tem)%nrec(2,sf_dyn(jf_tem)%nbb), wp )
768         ztintb =  1. - ztinta
769         IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN    ! Computes slopes (here avt is used as workspace)
770            uslp (:,:,:) = ztintb * uslpdta (:,:,:,1)  + ztinta * uslpdta (:,:,:,2) 
771            vslp (:,:,:) = ztintb * vslpdta (:,:,:,1)  + ztinta * vslpdta (:,:,:,2) 
772            wslpi(:,:,:) = ztintb * wslpidta(:,:,:,1)  + ztinta * wslpidta(:,:,:,2) 
773            wslpj(:,:,:) = ztintb * wslpjdta(:,:,:,1)  + ztinta * wslpjdta(:,:,:,2) 
774         ENDIF
775      ELSE
776         zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fnow(:,:,:) * tmask(:,:,:)   ! temperature
777         zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fnow(:,:,:) * tmask(:,:,:)   ! salinity
778         avt(:,:,:)        = sf_dyn(jf_avt)%fnow(:,:,:) * tmask(:,:,:)   ! vertical diffusive coef.
779         CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj, Kbb, Kmm )
780         !
781         IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN    ! Computes slopes (here avt is used as workspace)
782            uslp (:,:,:) = zuslp (:,:,:)
783            vslp (:,:,:) = zvslp (:,:,:)
784            wslpi(:,:,:) = zwslpi(:,:,:)
785            wslpj(:,:,:) = zwslpj(:,:,:)
786         ENDIF
787      ENDIF
788      !
789   END SUBROUTINE dta_dyn_slp
790
791
792   SUBROUTINE compute_slopes( kt, pts, puslp, pvslp, pwslpi, pwslpj, Kbb, Kmm )
793      !!---------------------------------------------------------------------
794      !!                    ***  ROUTINE dta_dyn_slp  ***
795      !!
796      !! ** Purpose :   Computation of slope
797      !!---------------------------------------------------------------------
798      INTEGER ,                              INTENT(in ) :: kt       ! time step
799      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in ) :: pts      ! temperature/salinity
800      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: puslp    ! zonal isopycnal slopes
801      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pvslp    ! meridional isopycnal slopes
802      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pwslpi   ! zonal diapycnal slopes
803      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pwslpj   ! meridional diapycnal slopes
804      INTEGER ,                              INTENT(in ) :: Kbb, Kmm ! ocean time level indices
805      !!---------------------------------------------------------------------
806      !
807      IF( l_ldfslp .AND. .NOT.lk_c1d ) THEN    ! Computes slopes (here avt is used as workspace)
808         CALL eos    ( pts, rhd, rhop, gdept_0(:,:,:) )
809         CALL eos_rab( pts, rab_n, Kmm )       ! now local thermal/haline expension ratio at T-points
810         CALL bn2    ( pts, rab_n, rn2, Kmm  ) ! now    Brunt-Vaisala
811
812      ! Partial steps: before Horizontal DErivative
813      IF( ln_zps  .AND. .NOT. ln_isfcav)                            &
814         &            CALL zps_hde    ( kt, Kmm, jpts, pts, gtsu, gtsv,  &  ! Partial steps: before horizontal gradient
815         &                                        rhd, gru , grv    )  ! of t, s, rd at the last ocean level
816      IF( ln_zps .AND.        ln_isfcav)                            &
817         &            CALL zps_hde_isf( kt, Kmm, jpts, pts, gtsu, gtsv, gtui, gtvi, &  ! Partial steps for top cell (ISF)
818         &                                        rhd, gru , grv , grui, grvi )  ! of t, s, rd at the first ocean level
819
820         rn2b(:,:,:) = rn2(:,:,:)                ! needed for zdfmxl
821         CALL zdf_mxl( kt, Kmm )                 ! mixed layer depth
822         CALL ldf_slp( kt, rhd, rn2, Kbb, Kmm )  ! slopes
823         puslp (:,:,:) = uslp (:,:,:)
824         pvslp (:,:,:) = vslp (:,:,:)
825         pwslpi(:,:,:) = wslpi(:,:,:)
826         pwslpj(:,:,:) = wslpj(:,:,:)
827     ELSE
828         puslp (:,:,:) = 0.            ! to avoid warning when compiling
829         pvslp (:,:,:) = 0.
830         pwslpi(:,:,:) = 0.
831         pwslpj(:,:,:) = 0.
832     ENDIF
833      !
834   END SUBROUTINE compute_slopes
835
836   !!======================================================================
837END MODULE dtadyn
Note: See TracBrowser for help on using the repository browser.