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/2013/dev_r3853_CNRS9_ConfSetting/NEMOGCM/NEMO/OFF_SRC – NEMO

source: branches/2013/dev_r3853_CNRS9_ConfSetting/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 3875

Last change on this file since 3875 was 3875, checked in by clevy, 11 years ago

Configuration Setting/Step? 1, see ticket:#1074

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