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/2012/dev_MERGE_2012/NEMOGCM/NEMO/OFF_SRC – NEMO

source: branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 3680

Last change on this file since 3680 was 3680, checked in by rblod, 11 years ago

First commit of the final branch for 2012 (future nemo_3_5), see ticket #1028

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