source: branches/2015/dev_r5204_CNRS_PISCES_dcy/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 5209

Last change on this file since 5209 was 5209, checked in by cetlod, 6 years ago

dev_r5204_CNRS_PISCES_dcy:updates for offline

  • Property svn:keywords set to Id
File size: 31.7 KB
Line 
1MODULE dtadyn
2   !!======================================================================
3   !!                       ***  MODULE  dtadyn  ***
4   !! Off-line : interpolation of the physical fields
5   !!======================================================================
6   !! History :   OPA  ! 1992-01 (M. Imbard) Original code
7   !!             8.0  ! 1998-04 (L.Bopp MA Foujols) slopes for isopyc.
8   !!              -   ! 1998-05 (L. Bopp) read output of coupled run
9   !!             8.2  ! 2001-01 (M. Levy et M. Benjelloul) add netcdf FORMAT
10   !!   NEMO      1.0  ! 2005-03 (O. Aumont and A. El Moussaoui) F90
11   !!              -   ! 2005-12 (C. Ethe) Adapted for DEGINT
12   !!             3.0  ! 2007-06 (C. Ethe) use of iom module
13   !!             3.3  ! 2010-11 (C. Ethe) Full reorganization of the off-line: phasing with the on-line
14   !!             3.4  ! 2011-05 (C. Ethe) Use of fldread
15   !!----------------------------------------------------------------------
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 zdf_oce         ! ocean vertical physics: variables
25   USE sbc_oce         ! surface module: variables
26   USE trc_oce         ! share ocean/biogeo variables
27   USE phycst          ! physical constants
28   USE trabbl          ! active tracer: bottom boundary layer
29   USE ldfslp          ! lateral diffusion: iso-neutral slopes
30   USE ldfeiv          ! eddy induced velocity coef.
31   USE ldftra_oce      ! ocean tracer   lateral physics
32   USE zdfmxl          ! vertical physics: mixed layer depth
33   USE eosbn2          ! equation of state - Brunt Vaisala frequency
34   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
35   USE zpshde          ! z-coord. with partial steps: horizontal derivatives
36   USE in_out_manager  ! I/O manager
37   USE iom             ! I/O library
38   USE lib_mpp         ! distributed memory computing library
39   USE prtctl          ! print control
40   USE fldread         ! read input fields
41   USE timing          ! Timing
42
43   IMPLICIT NONE
44   PRIVATE
45
46   PUBLIC   dta_dyn_init   ! called by opa.F90
47   PUBLIC   dta_dyn        ! called by step.F90
48
49   CHARACTER(len=100) ::   cn_dir       !: Root directory for location of ssr files
50   LOGICAL            ::   ln_dynwzv    !: vertical velocity read in a file (T) or computed from u/v (F)
51   LOGICAL            ::   ln_dynbbl    !: bbl coef read in a file (T) or computed (F)
52   LOGICAL            ::   ln_degrad    !: degradation option enabled or not
53   LOGICAL            ::   ln_dynrnf    !: read runoff data in file (T) or set to zero (F)
54
55   INTEGER  , PARAMETER ::   jpfld = 22     ! maximum number of fields to read
56   INTEGER  , SAVE      ::   jf_tem         ! index of temperature
57   INTEGER  , SAVE      ::   jf_sal         ! index of salinity
58   INTEGER  , SAVE      ::   jf_uwd         ! index of u-wind
59   INTEGER  , SAVE      ::   jf_vwd         ! index of v-wind
60   INTEGER  , SAVE      ::   jf_wwd         ! index of w-wind
61   INTEGER  , SAVE      ::   jf_avt         ! index of Kz
62   INTEGER  , SAVE      ::   jf_mld         ! index of mixed layer deptht
63   INTEGER  , SAVE      ::   jf_emp         ! index of water flux
64   INTEGER  , SAVE      ::   jf_qsr         ! index of solar radiation
65   INTEGER  , SAVE      ::   jf_qsm        ! index of daily mean solar radiation in case of dirunal cycle
66   INTEGER  , SAVE      ::   jf_wnd         ! index of wind speed
67   INTEGER  , SAVE      ::   jf_ice         ! index of sea ice cover
68   INTEGER  , SAVE      ::   jf_rnf         ! index of river runoff
69   INTEGER  , SAVE      ::   jf_ubl         ! index of u-bbl coef
70   INTEGER  , SAVE      ::   jf_vbl         ! index of v-bbl coef
71   INTEGER  , SAVE      ::   jf_ahu         ! index of u-diffusivity coef
72   INTEGER  , SAVE      ::   jf_ahv         ! index of v-diffusivity coef
73   INTEGER  , SAVE      ::   jf_ahw         ! index of w-diffusivity coef
74   INTEGER  , SAVE      ::   jf_eiu         ! index of u-eiv
75   INTEGER  , SAVE      ::   jf_eiv         ! index of v-eiv
76   INTEGER  , SAVE      ::   jf_eiw         ! index of w-eiv
77   INTEGER  , SAVE      ::   jf_fmf         ! index of downward salt flux
78
79   TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_dyn  ! structure of input fields (file informations, fields read)
80   !                                               !
81   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wdta       ! vertical velocity at 2 time step
82   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:  ) :: wnow       ! vertical velocity at 2 time step
83   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: uslpdta    ! zonal isopycnal slopes
84   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: vslpdta    ! meridional isopycnal slopes
85   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpidta   ! zonal diapycnal slopes
86   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpjdta   ! meridional diapycnal slopes
87   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:)   :: uslpnow    ! zonal isopycnal slopes
88   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:)   :: vslpnow    ! meridional isopycnal slopes
89   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:)   :: wslpinow   ! zonal diapycnal slopes
90   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:)   :: wslpjnow   ! meridional diapycnal slopes
91
92   INTEGER :: nrecprev_tem , nrecprev_uwd
93
94   !! * Substitutions
95#  include "domzgr_substitute.h90"
96#  include "vectopt_loop_substitute.h90"
97   !!----------------------------------------------------------------------
98   !! NEMO/OFF 3.3 , NEMO Consortium (2010)
99   !! $Id$
100   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
101   !!----------------------------------------------------------------------
102CONTAINS
103
104   SUBROUTINE dta_dyn( kt )
105      !!----------------------------------------------------------------------
106      !!                  ***  ROUTINE dta_dyn  ***
107      !!
108      !! ** Purpose :  Prepares dynamics and physics fields from a NEMO run
109      !!               for an off-line simulation of passive tracers
110      !!
111      !! ** Method : calculates the position of data
112      !!             - computes slopes if needed
113      !!             - interpolates data if needed
114      !!----------------------------------------------------------------------
115      !
116      USE oce, ONLY:  zts    => tsa 
117      USE oce, ONLY:  zuslp  => ua   , zvslp  => va
118      USE oce, ONLY:  zwslpi => rotb , zwslpj => rotn
119      USE oce, ONLY:  zu     => ub   , zv     => vb,  zw => hdivb
120      !
121      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
122      !
123      INTEGER  ::   ji, jj     ! dummy loop indices
124      INTEGER  ::   isecsbc    ! number of seconds between Jan. 1st 00h of nit000 year and the middle of time step
125      REAL(wp) ::   ztinta     ! ratio applied to after  records when doing time interpolation
126      REAL(wp) ::   ztintb     ! ratio applied to before records when doing time interpolation
127      INTEGER  ::   iswap_tem, iswap_uwd    !
128      !!----------------------------------------------------------------------
129     
130      !
131      IF( nn_timing == 1 )  CALL timing_start( 'dta_dyn')
132      !
133      isecsbc = nsec_year + nsec1jan000 
134      !
135      IF( kt == nit000 ) THEN
136         nrecprev_tem = 0
137         nrecprev_uwd = 0
138         !
139         CALL fld_read( kt, 1, sf_dyn )      !==   read data at kt time step   ==!
140         !
141         IF( lk_ldfslp .AND. .NOT.lk_c1d .AND. sf_dyn(jf_tem)%ln_tint ) THEN    ! Computes slopes (here avt is used as workspace)                       
142            zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,1) * tmask(:,:,:)   ! temperature
143            zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,1) * tmask(:,:,:)   ! salinity
144            avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,1) * tmask(:,:,:)   ! vertical diffusive coef.
145            CALL dta_dyn_slp( kt, zts, zuslp, zvslp, zwslpi, zwslpj )
146            uslpdta (:,:,:,1) = zuslp (:,:,:) 
147            vslpdta (:,:,:,1) = zvslp (:,:,:) 
148            wslpidta(:,:,:,1) = zwslpi(:,:,:) 
149            wslpjdta(:,:,:,1) = zwslpj(:,:,:) 
150         ENDIF
151         IF( ln_dynwzv .AND. sf_dyn(jf_uwd)%ln_tint )  THEN    ! compute vertical velocity from u/v
152            zu(:,:,:) = sf_dyn(jf_uwd)%fdta(:,:,:,1)
153            zv(:,:,:) = sf_dyn(jf_vwd)%fdta(:,:,:,1)
154            CALL dta_dyn_wzv( zu, zv, zw )
155            wdta(:,:,:,1) = zw(:,:,:) * tmask(:,:,:)
156         ENDIF
157      ELSE
158         nrecprev_tem = sf_dyn(jf_tem)%nrec_a(2)
159         nrecprev_uwd = sf_dyn(jf_uwd)%nrec_a(2)
160         !
161         CALL fld_read( kt, 1, sf_dyn )      !==   read data at kt time step   ==!
162         !
163      ENDIF
164      !
165      IF( lk_ldfslp .AND. .NOT.lk_c1d ) THEN    ! Computes slopes (here avt is used as workspace)                       
166         iswap_tem = 0
167         IF(  kt /= nit000 .AND. ( sf_dyn(jf_tem)%nrec_a(2) - nrecprev_tem ) /= 0 )  iswap_tem = 1
168         IF( ( isecsbc > sf_dyn(jf_tem)%nrec_b(2) .AND. iswap_tem == 1 ) .OR. kt == nit000 )  THEN    ! read/update the after data
169            IF(lwp) WRITE(numout,*) ' Compute new slopes at kt = ', kt
170            IF( sf_dyn(jf_tem)%ln_tint ) THEN                 ! time interpolation of data
171               IF( kt /= nit000 ) THEN
172                  uslpdta (:,:,:,1) =  uslpdta (:,:,:,2)         ! swap the data
173                  vslpdta (:,:,:,1) =  vslpdta (:,:,:,2) 
174                  wslpidta(:,:,:,1) =  wslpidta(:,:,:,2) 
175                  wslpjdta(:,:,:,1) =  wslpjdta(:,:,:,2) 
176               ENDIF
177               !
178               zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,2) * tmask(:,:,:)   ! temperature
179               zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,2) * tmask(:,:,:)   ! salinity
180               avt(:,:,:)        = sf_dyn(jf_avt)%fdta(:,:,:,2) * tmask(:,:,:)   ! vertical diffusive coef.
181               CALL dta_dyn_slp( kt, zts, zuslp, zvslp, zwslpi, zwslpj )
182               !
183               uslpdta (:,:,:,2) = zuslp (:,:,:) 
184               vslpdta (:,:,:,2) = zvslp (:,:,:) 
185               wslpidta(:,:,:,2) = zwslpi(:,:,:) 
186               wslpjdta(:,:,:,2) = zwslpj(:,:,:) 
187            ELSE
188               zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fnow(:,:,:) * tmask(:,:,:)
189               zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fnow(:,:,:) * tmask(:,:,:)
190               avt(:,:,:)        = sf_dyn(jf_avt)%fnow(:,:,:) * tmask(:,:,:)
191               CALL dta_dyn_slp( kt, zts, zuslp, zvslp, zwslpi, zwslpj )
192               uslpnow (:,:,:)   = zuslp (:,:,:) 
193               vslpnow (:,:,:)   = zvslp (:,:,:) 
194               wslpinow(:,:,:)   = zwslpi(:,:,:) 
195               wslpjnow(:,:,:)   = zwslpj(:,:,:) 
196            ENDIF
197         ENDIF
198         IF( sf_dyn(jf_tem)%ln_tint )  THEN
199            ztinta =  REAL( isecsbc - sf_dyn(jf_tem)%nrec_b(2), wp )  &
200               &    / REAL( sf_dyn(jf_tem)%nrec_a(2) - sf_dyn(jf_tem)%nrec_b(2), wp )
201            ztintb =  1. - ztinta
202            uslp (:,:,:) = ztintb * uslpdta (:,:,:,1)  + ztinta * uslpdta (:,:,:,2) 
203            vslp (:,:,:) = ztintb * vslpdta (:,:,:,1)  + ztinta * vslpdta (:,:,:,2) 
204            wslpi(:,:,:) = ztintb * wslpidta(:,:,:,1)  + ztinta * wslpidta(:,:,:,2) 
205            wslpj(:,:,:) = ztintb * wslpjdta(:,:,:,1)  + ztinta * wslpjdta(:,:,:,2) 
206         ELSE
207            uslp (:,:,:) = uslpnow (:,:,:)
208            vslp (:,:,:) = vslpnow (:,:,:)
209            wslpi(:,:,:) = wslpinow(:,:,:)
210            wslpj(:,:,:) = wslpjnow(:,:,:)
211         ENDIF
212      ENDIF
213      !
214      IF( ln_dynwzv )  THEN    ! compute vertical velocity from u/v
215         iswap_uwd = 0
216         IF(  kt /= nit000 .AND. ( sf_dyn(jf_uwd)%nrec_a(2) - nrecprev_uwd ) /= 0 )  iswap_uwd = 1
217         IF( ( isecsbc > sf_dyn(jf_uwd)%nrec_b(2) .AND. iswap_uwd == 1 ) .OR. kt == nit000 )  THEN    ! read/update the after data
218            IF(lwp) WRITE(numout,*) ' Compute new vertical velocity at kt = ', kt
219            IF(lwp) WRITE(numout,*)
220            IF( sf_dyn(jf_uwd)%ln_tint ) THEN                 ! time interpolation of data
221               IF( kt /= nit000 )  THEN
222                  wdta(:,:,:,1) =  wdta(:,:,:,2)     ! swap the data for initialisation
223               ENDIF
224               zu(:,:,:) = sf_dyn(jf_uwd)%fdta(:,:,:,2)
225               zv(:,:,:) = sf_dyn(jf_vwd)%fdta(:,:,:,2)
226               CALL dta_dyn_wzv( zu, zv, zw )
227               wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:)
228            ELSE
229               zu(:,:,:) = sf_dyn(jf_uwd)%fnow(:,:,:) 
230               zv(:,:,:) = sf_dyn(jf_vwd)%fnow(:,:,:)
231               CALL dta_dyn_wzv( zu, zv, zw )
232               wnow(:,:,:)  = zw(:,:,:) * tmask(:,:,:)
233            ENDIF
234         ENDIF
235         IF( sf_dyn(jf_uwd)%ln_tint )  THEN
236            ztinta =  REAL( isecsbc - sf_dyn(jf_uwd)%nrec_b(2), wp )  &
237               &    / REAL( sf_dyn(jf_uwd)%nrec_a(2) - sf_dyn(jf_uwd)%nrec_b(2), wp )
238            ztintb =  1. - ztinta
239            wn(:,:,:) = ztintb * wdta(:,:,:,1)  + ztinta * wdta(:,:,:,2) 
240         ELSE
241            wn(:,:,:) = wnow(:,:,:)
242         ENDIF
243      ENDIF
244      !
245      tsn(:,:,:,jp_tem) = sf_dyn(jf_tem)%fnow(:,:,:)  * tmask(:,:,:)    ! temperature
246      tsn(:,:,:,jp_sal) = sf_dyn(jf_sal)%fnow(:,:,:)  * tmask(:,:,:)    ! salinity
247      !
248      !
249      CALL eos    ( tsn, rhd, rhop, gdept_0(:,:,:) ) ! In any case, we need rhop
250      CALL eos_rab( tsn, rab_n )       ! now    local thermal/haline expension ratio at T-points
251      CALL bn2    ( tsn, rab_n, rn2 ) ! before Brunt-Vaisala frequency need for zdfmxl
252
253      rn2b(:,:,:) = rn2(:,:,:)         ! need for zdfmxl
254      CALL zdf_mxl( kt )                                                   ! In any case, we need mxl
255      !
256      avt(:,:,:)       = sf_dyn(jf_avt)%fnow(:,:,:)  * tmask(:,:,:)    ! vertical diffusive coefficient
257      un (:,:,:)       = sf_dyn(jf_uwd)%fnow(:,:,:)  * umask(:,:,:)    ! u-velocity
258      vn (:,:,:)       = sf_dyn(jf_vwd)%fnow(:,:,:)  * vmask(:,:,:)    ! v-velocity
259      IF( .NOT.ln_dynwzv ) &                                          ! w-velocity read in file
260         wn (:,:,:)    = sf_dyn(jf_wwd)%fnow(:,:,:) * tmask(:,:,:)   
261      hmld(:,:)        = sf_dyn(jf_mld)%fnow(:,:,1) * tmask(:,:,1)    ! mixed layer depht
262      wndm(:,:)        = sf_dyn(jf_wnd)%fnow(:,:,1) * tmask(:,:,1)    ! wind speed - needed for gas exchange
263      emp (:,:)        = sf_dyn(jf_emp)%fnow(:,:,1) * tmask(:,:,1)    ! E-P
264      fmmflx(:,:)      = sf_dyn(jf_fmf)%fnow(:,:,1) * tmask(:,:,1)    ! downward salt flux (v3.5+)
265      fr_i(:,:)        = sf_dyn(jf_ice)%fnow(:,:,1) * tmask(:,:,1)    ! Sea-ice fraction
266      qsr (:,:)        = sf_dyn(jf_qsr)%fnow(:,:,1) * tmask(:,:,1)    ! solar radiation
267      IF( ln_dm2dc ) THEN  ! diurnal cycle
268         qsr_mean(:,:) = sf_dyn(jf_qsm)%fnow(:,:,1) * tmask(:,:,1)    ! daily mean solar radiation
269      ELSE
270         qsr_mean(:,:) = qsr(:,:)   ! solar radiation
271      ENDIF
272      IF( ln_dynrnf ) &
273      rnf (:,:)        = sf_dyn(jf_rnf)%fnow(:,:,1) * tmask(:,:,1)    ! river runoffs
274
275      !                                                      ! bbl diffusive coef
276#if defined key_trabbl && ! defined key_c1d
277      IF( ln_dynbbl ) THEN                                        ! read in a file
278         ahu_bbl(:,:)  = sf_dyn(jf_ubl)%fnow(:,:,1) * umask(:,:,1)
279         ahv_bbl(:,:)  = sf_dyn(jf_vbl)%fnow(:,:,1) * vmask(:,:,1)
280      ELSE                                                        ! Compute bbl coefficients if needed
281         tsb(:,:,:,:) = tsn(:,:,:,:)
282         CALL bbl( kt, nit000, 'TRC')
283      END IF
284#endif
285#if ( ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv ) && ! defined key_c1d
286      aeiw(:,:)        = sf_dyn(jf_eiw)%fnow(:,:,1) * tmask(:,:,1)    ! w-eiv
287      !                                                           ! Computes the horizontal values from the vertical value
288      DO jj = 2, jpjm1
289         DO ji = fs_2, fs_jpim1   ! vector opt.
290            aeiu(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji+1,jj  ) )  ! Average the diffusive coefficient at u- v- points
291            aeiv(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji  ,jj+1) )  ! at u- v- points
292         END DO
293      END DO
294      CALL lbc_lnk( aeiu, 'U', 1. )   ;   CALL lbc_lnk( aeiv, 'V', 1. )    ! lateral boundary condition
295#endif
296     
297#if defined key_degrad && ! defined key_c1d
298      !                                          ! degrad option : diffusive and eiv coef are 3D
299      ahtu(:,:,:) = sf_dyn(jf_ahu)%fnow(:,:,:) * umask(:,:,:)
300      ahtv(:,:,:) = sf_dyn(jf_ahv)%fnow(:,:,:) * vmask(:,:,:)
301      ahtw(:,:,:) = sf_dyn(jf_ahw)%fnow(:,:,:) * tmask(:,:,:)
302#  if defined key_traldf_eiv 
303      aeiu(:,:,:) = sf_dyn(jf_eiu)%fnow(:,:,:) * umask(:,:,:)
304      aeiv(:,:,:) = sf_dyn(jf_eiv)%fnow(:,:,:) * vmask(:,:,:)
305      aeiw(:,:,:) = sf_dyn(jf_eiw)%fnow(:,:,:) * tmask(:,:,:)
306#  endif
307#endif
308      !
309      IF(ln_ctl) THEN                  ! print control
310         CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' tn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
311         CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_sal), clinfo1=' sn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
312         CALL prt_ctl(tab3d_1=un               , clinfo1=' un      - : ', mask1=umask, ovlap=1, kdim=jpk   )
313         CALL prt_ctl(tab3d_1=vn               , clinfo1=' vn      - : ', mask1=vmask, ovlap=1, kdim=jpk   )
314         CALL prt_ctl(tab3d_1=wn               , clinfo1=' wn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
315         CALL prt_ctl(tab3d_1=avt              , clinfo1=' kz      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
316         CALL prt_ctl(tab2d_1=fr_i             , clinfo1=' fr_i    - : ', mask1=tmask, ovlap=1 )
317         CALL prt_ctl(tab2d_1=hmld             , clinfo1=' hmld    - : ', mask1=tmask, ovlap=1 )
318         CALL prt_ctl(tab2d_1=fmmflx           , clinfo1=' fmmflx  - : ', mask1=tmask, ovlap=1 )
319         CALL prt_ctl(tab2d_1=emp              , clinfo1=' emp     - : ', mask1=tmask, ovlap=1 )
320         CALL prt_ctl(tab2d_1=wndm             , clinfo1=' wspd    - : ', mask1=tmask, ovlap=1 )
321         CALL prt_ctl(tab2d_1=qsr              , clinfo1=' qsr     - : ', mask1=tmask, ovlap=1 )
322      ENDIF
323      !
324      IF( nn_timing == 1 )  CALL timing_stop( 'dta_dyn')
325      !
326   END SUBROUTINE dta_dyn
327
328
329   SUBROUTINE dta_dyn_init
330      !!----------------------------------------------------------------------
331      !!                  ***  ROUTINE dta_dyn_init  ***
332      !!
333      !! ** Purpose :   Initialisation of the dynamical data     
334      !! ** Method  : - read the data namdta_dyn namelist
335      !!
336      !! ** Action  : - read parameters
337      !!----------------------------------------------------------------------
338      INTEGER  :: ierr, ierr0, ierr1, ierr2, ierr3   ! return error code
339      INTEGER  :: ifpr                               ! dummy loop indice
340      INTEGER  :: jfld                               ! dummy loop arguments
341      INTEGER  :: inum, idv, idimv                   ! local integer
342      INTEGER  :: ios                                ! Local integer output status for namelist read
343      !!
344      CHARACTER(len=100)            ::  cn_dir   !   Root directory for location of core files
345      TYPE(FLD_N), DIMENSION(jpfld) ::  slf_d    ! array of namelist informations on the fields to read
346      TYPE(FLD_N) :: sn_tem, sn_sal, sn_mld, sn_emp, sn_ice, sn_qsr, sn_qsm, sn_wnd, sn_rnf  ! informations about the fields to be read
347      TYPE(FLD_N) :: sn_uwd, sn_vwd, sn_wwd, sn_avt, sn_ubl, sn_vbl          !   "                                 "
348      TYPE(FLD_N) :: sn_ahu, sn_ahv, sn_ahw, sn_eiu, sn_eiv, sn_eiw, sn_fmf  !   "                                 "
349      !!----------------------------------------------------------------------
350      !
351      NAMELIST/namdta_dyn/cn_dir, ln_dynwzv, ln_dynbbl, ln_degrad, ln_dynrnf,    &
352         &                sn_tem, sn_sal, sn_mld, sn_emp, sn_ice, sn_qsr, sn_qsm, sn_wnd, sn_rnf,  &
353         &                sn_uwd, sn_vwd, sn_wwd, sn_avt, sn_ubl, sn_vbl,          &
354         &                sn_ahu, sn_ahv, sn_ahw, sn_eiu, sn_eiv, sn_eiw, sn_fmf
355      !
356      REWIND( numnam_ref )              ! Namelist namdta_dyn in reference namelist : Offline: init. of dynamical data
357      READ  ( numnam_ref, namdta_dyn, IOSTAT = ios, ERR = 901)
358901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdta_dyn in reference namelist', lwp )
359
360      REWIND( numnam_cfg )              ! Namelist namdta_dyn in configuration namelist : Offline: init. of dynamical data
361      READ  ( numnam_cfg, namdta_dyn, IOSTAT = ios, ERR = 902 )
362902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdta_dyn in configuration namelist', lwp )
363      IF(lwm) WRITE ( numond, namdta_dyn )
364      !                                         ! store namelist information in an array
365      !                                         ! Control print
366      IF(lwp) THEN
367         WRITE(numout,*)
368         WRITE(numout,*) 'dta_dyn : offline dynamics '
369         WRITE(numout,*) '~~~~~~~ '
370         WRITE(numout,*) '   Namelist namdta_dyn'
371         WRITE(numout,*) '      vertical velocity read from file (T) or computed (F) ln_dynwzv  = ', ln_dynwzv
372         WRITE(numout,*) '      bbl coef read from file (T) or computed (F)          ln_dynbbl  = ', ln_dynbbl
373         WRITE(numout,*) '      degradation option enabled (T) or not (F)            ln_degrad  = ', ln_degrad
374         WRITE(numout,*) '      river runoff option enabled (T) or not (F)           ln_dynrnf  = ', ln_dynrnf
375         WRITE(numout,*)
376      ENDIF
377      !
378      IF( ln_degrad .AND. .NOT.lk_degrad ) THEN
379         CALL ctl_warn( 'dta_dyn_init: degradation option requires key_degrad activated ; force ln_degrad to false' )
380         ln_degrad = .FALSE.
381      ENDIF
382      IF( ln_dynbbl .AND. ( .NOT.lk_trabbl .OR. lk_c1d ) ) THEN
383         CALL ctl_warn( 'dta_dyn_init: bbl option requires key_trabbl activated ; force ln_dynbbl to false' )
384         ln_dynbbl = .FALSE.
385      ENDIF
386
387      jf_tem = 1   ;   jf_sal = 2   ;  jf_mld = 3   ;  jf_emp = 4   ;   jf_fmf  = 5   ;  jf_ice = 6   ;   jf_qsr = 7
388      jf_wnd = 8   ;   jf_uwd = 9   ;  jf_vwd = 10  ;  jf_wwd = 11  ;   jf_avt  = 12  ;  jfld  = jf_avt
389      !
390      slf_d(jf_tem) = sn_tem   ;   slf_d(jf_sal)  = sn_sal   ;   slf_d(jf_mld) = sn_mld
391      slf_d(jf_emp) = sn_emp   ;   slf_d(jf_fmf ) = sn_fmf   ;   slf_d(jf_ice) = sn_ice 
392      slf_d(jf_qsr) = sn_qsr   ;   slf_d(jf_wnd)  = sn_wnd   ;   slf_d(jf_avt) = sn_avt 
393      slf_d(jf_uwd) = sn_uwd   ;   slf_d(jf_vwd)  = sn_vwd   ;   slf_d(jf_wwd) = sn_wwd
394
395      !
396      IF( ln_dynrnf ) THEN
397                jf_rnf = jfld + 1  ;  jfld  = jf_rnf
398         slf_d(jf_rnf) = sn_rnf
399      ELSE
400         rnf (:,:) = 0._wp
401      ENDIF
402
403      IF( ln_dm2dc ) THEN
404               jf_qsm  = jfld + 1  ;  jfld  = jf_qsm
405         slf_d(jf_qsm) = sn_qsm
406      ENDIF
407
408      !
409      IF( .NOT.ln_degrad ) THEN     ! no degrad option
410         IF( lk_traldf_eiv .AND. ln_dynbbl ) THEN        ! eiv & bbl
411                 jf_ubl  = jfld + 1 ;        jf_vbl  = jfld + 2 ;        jf_eiw  = jfld + 3   ;   jfld = jf_eiw
412           slf_d(jf_ubl) = sn_ubl  ;   slf_d(jf_vbl) = sn_vbl  ;   slf_d(jf_eiw) = sn_eiw
413         ENDIF
414         IF( .NOT.lk_traldf_eiv .AND. ln_dynbbl ) THEN   ! no eiv & bbl
415                 jf_ubl  = jfld + 1 ;        jf_vbl  = jfld + 2 ;  jfld = jf_vbl
416           slf_d(jf_ubl) = sn_ubl  ;   slf_d(jf_vbl) = sn_vbl
417         ENDIF
418         IF( lk_traldf_eiv .AND. .NOT.ln_dynbbl ) THEN   ! eiv & no bbl
419           jf_eiw = jfld + 1 ; jfld = jf_eiw ; slf_d(jf_eiw) = sn_eiw
420         ENDIF
421      ELSE
422              jf_ahu  = jfld + 1 ;        jf_ahv  = jfld + 2 ;        jf_ahw  = jfld + 3  ;  jfld = jf_ahw
423        slf_d(jf_ahu) = sn_ahu  ;   slf_d(jf_ahv) = sn_ahv  ;   slf_d(jf_ahw) = sn_ahw
424        IF( lk_traldf_eiv .AND. ln_dynbbl ) THEN         ! eiv & bbl
425                 jf_ubl  = jfld + 1 ;        jf_vbl  = jfld + 2 ;
426           slf_d(jf_ubl) = sn_ubl  ;   slf_d(jf_vbl) = sn_vbl
427                 jf_eiu  = jfld + 3 ;        jf_eiv  = jfld + 4 ;    jf_eiw  = jfld + 5   ;  jfld = jf_eiw 
428           slf_d(jf_eiu) = sn_eiu  ;   slf_d(jf_eiv) = sn_eiv  ;    slf_d(jf_eiw) = sn_eiw
429        ENDIF
430        IF( .NOT.lk_traldf_eiv .AND. ln_dynbbl ) THEN    ! no eiv & bbl
431                 jf_ubl  = jfld + 1 ;        jf_vbl  = jfld + 2 ;  jfld = jf_vbl
432           slf_d(jf_ubl) = sn_ubl  ;   slf_d(jf_vbl) = sn_vbl
433        ENDIF
434        IF( lk_traldf_eiv .AND. .NOT.ln_dynbbl ) THEN    ! eiv & no bbl
435                 jf_eiu  = jfld + 1 ;         jf_eiv  = jfld + 2 ;    jf_eiw  = jfld + 3   ; jfld = jf_eiw 
436           slf_d(jf_eiu) = sn_eiu  ;   slf_d(jf_eiv) = sn_eiv  ;   slf_d(jf_eiw) = sn_eiw
437        ENDIF
438      ENDIF
439 
440      ALLOCATE( sf_dyn(jfld), STAT=ierr )         ! set sf structure
441      IF( ierr > 0 ) THEN
442         CALL ctl_stop( 'dta_dyn: unable to allocate sf structure' )   ;   RETURN
443      ENDIF
444      ! Open file for each variable to get his number of dimension
445      DO ifpr = 1, jfld
446         CALL iom_open( TRIM( cn_dir )//TRIM( slf_d(ifpr)%clname ), inum )
447         idv   = iom_varid( inum , slf_d(ifpr)%clvar )  ! id of the variable sdjf%clvar
448         idimv = iom_file ( inum )%ndims(idv)             ! number of dimension for variable sdjf%clvar
449         IF( inum /= 0 )   CALL iom_close( inum )       ! close file if already open
450         IF( idimv == 3 ) THEN    ! 2D variable
451                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,1)    , STAT=ierr0 )
452            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,1,2)  , STAT=ierr1 )
453         ELSE                     ! 3D variable
454                                      ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,jpk)  , STAT=ierr0 )
455            IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,jpk,2), STAT=ierr1 )
456         ENDIF
457         IF( ierr0 + ierr1 > 0 ) THEN
458            CALL ctl_stop( 'dta_dyn_init : unable to allocate sf_dyn array structure' )   ;   RETURN
459         ENDIF
460      END DO
461      !                                         ! fill sf with slf_i and control print
462      CALL fld_fill( sf_dyn, slf_d, cn_dir, 'dta_dyn_init', 'Data in file', 'namdta_dyn' )
463      !
464      IF( lk_ldfslp .AND. .NOT.lk_c1d ) THEN                  ! slopes
465         IF( sf_dyn(jf_tem)%ln_tint ) THEN      ! time interpolation
466            ALLOCATE( uslpdta (jpi,jpj,jpk,2), vslpdta (jpi,jpj,jpk,2),    &
467            &         wslpidta(jpi,jpj,jpk,2), wslpjdta(jpi,jpj,jpk,2), STAT=ierr2 )
468         ELSE
469            ALLOCATE( uslpnow (jpi,jpj,jpk)  , vslpnow (jpi,jpj,jpk)  ,    &
470            &         wslpinow(jpi,jpj,jpk)  , wslpjnow(jpi,jpj,jpk)  , STAT=ierr2 )
471         ENDIF
472         IF( ierr2 > 0 ) THEN
473            CALL ctl_stop( 'dta_dyn_init : unable to allocate slope arrays' )   ;   RETURN
474         ENDIF
475      ENDIF
476      IF( ln_dynwzv ) THEN                  ! slopes
477         IF( sf_dyn(jf_uwd)%ln_tint ) THEN      ! time interpolation
478            ALLOCATE( wdta(jpi,jpj,jpk,2), STAT=ierr3 )
479         ELSE
480            ALLOCATE( wnow(jpi,jpj,jpk)  , STAT=ierr3 )
481         ENDIF
482         IF( ierr3 > 0 ) THEN
483            CALL ctl_stop( 'dta_dyn_init : unable to allocate wdta arrays' )   ;   RETURN
484         ENDIF
485      ENDIF
486      !
487      CALL dta_dyn( nit000 )
488      !
489   END SUBROUTINE dta_dyn_init
490
491   SUBROUTINE dta_dyn_wzv( pu, pv, pw )
492      !!----------------------------------------------------------------------
493      !!                    ***  ROUTINE wzv  ***
494      !!
495      !! ** Purpose :   Compute the now vertical velocity after the array swap
496      !!
497      !! ** Method  : - compute the now divergence given by :
498      !!         * z-coordinate ONLY !!!!
499      !!         hdiv = 1/(e1t*e2t) [ di(e2u  u) + dj(e1v  v) ]
500      !!     - Using the incompressibility hypothesis, the vertical
501      !!      velocity is computed by integrating the horizontal divergence
502      !!      from the bottom to the surface.
503      !!        The boundary conditions are w=0 at the bottom (no flux).
504      !!----------------------------------------------------------------------
505      USE oce, ONLY:  zhdiv => hdivn
506      !
507      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) :: pu, pv    !:  horizontal velocities
508      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(  out) :: pw        !:  vertical velocity
509      !!
510      INTEGER  ::  ji, jj, jk
511      REAL(wp) ::  zu, zu1, zv, zv1, zet
512      !!----------------------------------------------------------------------
513      !
514      ! Computation of vertical velocity using horizontal divergence
515      zhdiv(:,:,:) = 0._wp
516      DO jk = 1, jpkm1
517         DO jj = 2, jpjm1
518            DO ji = fs_2, fs_jpim1   ! vector opt.
519               zu  = pu(ji  ,jj  ,jk) * umask(ji  ,jj  ,jk) * e2u(ji  ,jj  ) * fse3u(ji  ,jj  ,jk)
520               zu1 = pu(ji-1,jj  ,jk) * umask(ji-1,jj  ,jk) * e2u(ji-1,jj  ) * fse3u(ji-1,jj  ,jk)
521               zv  = pv(ji  ,jj  ,jk) * vmask(ji  ,jj  ,jk) * e1v(ji  ,jj  ) * fse3v(ji  ,jj  ,jk)
522               zv1 = pv(ji  ,jj-1,jk) * vmask(ji  ,jj-1,jk) * e1v(ji  ,jj-1) * fse3v(ji  ,jj-1,jk)
523               zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
524               zhdiv(ji,jj,jk) = ( zu - zu1 + zv - zv1 ) * zet 
525            END DO
526         END DO
527      END DO
528      CALL lbc_lnk( zhdiv, 'T', 1. )      ! Lateral boundary conditions on zhdiv
529      !
530      ! computation of vertical velocity from the bottom
531      pw(:,:,jpk) = 0._wp
532      DO jk = jpkm1, 1, -1
533         pw(:,:,jk) = pw(:,:,jk+1) - fse3t(:,:,jk) * zhdiv(:,:,jk)
534      END DO
535      !
536   END SUBROUTINE dta_dyn_wzv
537
538   SUBROUTINE dta_dyn_slp( kt, pts, puslp, pvslp, pwslpi, pwslpj )
539      !!---------------------------------------------------------------------
540      !!                    ***  ROUTINE dta_dyn_slp  ***
541      !!
542      !! ** Purpose : Computation of slope
543      !!
544      !!---------------------------------------------------------------------
545      INTEGER ,                              INTENT(in ) :: kt       ! time step
546      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in ) :: pts      ! temperature/salinity
547      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: puslp    ! zonal isopycnal slopes
548      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pvslp    ! meridional isopycnal slopes
549      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pwslpi   ! zonal diapycnal slopes
550      REAL(wp), DIMENSION(jpi,jpj,jpk)     , INTENT(out) :: pwslpj   ! meridional diapycnal slopes
551      !!---------------------------------------------------------------------
552#if defined key_ldfslp && ! defined key_c1d
553      CALL eos    ( pts, rhd, rhop, gdept_0(:,:,:) )
554      CALL eos_rab( pts, rab_n )       ! now local thermal/haline expension ratio at T-points
555      CALL bn2    ( pts, rab_n, rn2  ) ! now    Brunt-Vaisala
556
557      ! Partial steps: before Horizontal DErivative
558      IF( ln_zps  .AND. .NOT. ln_isfcav)                            &
559         &            CALL zps_hde    ( kt, jpts, pts, gtsu, gtsv,  &  ! Partial steps: before horizontal gradient
560         &                                        rhd, gru , grv    )  ! of t, s, rd at the last ocean level
561      IF( ln_zps .AND.        ln_isfcav)                            &
562         &            CALL zps_hde_isf( kt, jpts, pts, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
563         &                                        rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
564         &                                 gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the first ocean level
565
566      rn2b(:,:,:) = rn2(:,:,:)         ! need for zdfmxl
567      CALL zdf_mxl( kt )            ! mixed layer depth
568      CALL ldf_slp( kt, rhd, rn2 )  ! slopes
569      puslp (:,:,:) = uslp (:,:,:) 
570      pvslp (:,:,:) = vslp (:,:,:) 
571      pwslpi(:,:,:) = wslpi(:,:,:) 
572      pwslpj(:,:,:) = wslpj(:,:,:) 
573#else
574      puslp (:,:,:) = 0.            ! to avoid warning when compiling
575      pvslp (:,:,:) = 0.
576      pwslpi(:,:,:) = 0.
577      pwslpj(:,:,:) = 0.
578#endif
579      !
580   END SUBROUTINE dta_dyn_slp
581   !!======================================================================
582END MODULE dtadyn
Note: See TracBrowser for help on using the repository browser.