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

source: NEMO/branches/UKMO/dev_r9885_GO6_mixing/src/OFF/dtadyn.F90 @ 9895

Last change on this file since 9895 was 9895, checked in by davestorkey, 6 years ago

UKMO/dev_r9885_GO6_mixing branch: clear SVN keywords.

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