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/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OFF_SRC – NEMO

source: branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 3119

Last change on this file since 3119 was 3119, checked in by cetlod, 12 years ago

branch dev_NEMO_MERGE_2011/NEMOGCM: minor bug corrections for offline configuration

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