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

Last change on this file was 15090, checked in by cetlod, 3 years ago

trunk : minor changes in PISCES to make it work in debug mode when nn_hls=2 ( still not working when enabling tiles )

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