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

source: branches/2011/dev_r2787_LOCEAN_offline_fldread/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 2821

Last change on this file since 2821 was 2821, checked in by cetlod, 13 years ago

Use of fldread for offline, see ticket #853

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