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 branches/2015/dev_r6946_Offline_vvl/NEMOGCM/NEMO/OFF_SRC – NEMO

source: branches/2015/dev_r6946_Offline_vvl/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 6952

Last change on this file since 6952 was 6952, checked in by cetlod, 8 years ago

Offline vvl: 1st implementation, see ticket #1775

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