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

source: trunk/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 2559

Last change on this file since 2559 was 2559, checked in by rblod, 13 years ago

Correction in dtadyn when not using key_trabbl

  • Property svn:keywords set to Id
File size: 34.1 KB
Line 
1MODULE dtadyn
2   !!======================================================================
3   !!                       ***  MODULE  dtadyn  ***
4   !! Off-line : interpolation of the physical fields
5   !!======================================================================
6   !! History :   OPA  ! 1992-01 (M. Imbard) Original code
7   !!             8.0  ! 1998-04 (L.Bopp MA Foujols) slopes for isopyc.
8   !!              -   ! 1998-05 (L. Bopp) read output of coupled run
9   !!             8.2  ! 2001-01 (M. Levy et M. Benjelloul) add netcdf FORMAT
10   !!   NEMO      1.0  ! 2005-03 (O. Aumont and A. El Moussaoui) F90
11   !!              -   ! 2005-12 (C. Ethe) Adapted for DEGINT
12   !!             3.0  ! 2007-06 (C. Ethe) use of iom module
13   !!              -   ! 2007-09  (C. Ethe)  add swap_dyn_data
14   !!             3.3  ! 2010-11 (C. Ethe) Full reorganization of the off-line: phasing with the on-line
15   !!----------------------------------------------------------------------
16
17   !!----------------------------------------------------------------------
18   !!   dta_dyn_init : initialization, namelist read, and parameters 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 phycst          ! physical constants
27   USE trabbl          ! active tracer: bottom boundary layer
28   USE ldfslp          ! lateral diffusion: iso-neutral slopes
29   USE ldfeiv          ! eddy induced velocity coef.
30   USE ldftra_oce      ! ocean tracer   lateral physics
31   USE zdfmxl          ! vertical physics: mixed layer depth
32   USE eosbn2          ! equation of state - Brunt Vaisala frequency
33   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
34   USE zpshde          ! z-coord. with partial steps: horizontal derivatives
35   USE in_out_manager  ! I/O manager
36   USE iom             ! I/O library
37   USE lib_mpp         ! distributed memory computing library
38   USE prtctl          !  print control
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   dta_dyn_init   ! called by opa.F90
44   PUBLIC   dta_dyn        ! called by step.F90
45
46   LOGICAL, PUBLIC ::   lperdyn = .TRUE.   !: boolean for periodic fields or not
47   LOGICAL, PUBLIC ::   lfirdyn = .TRUE.   !: boolean for the first call or not
48
49   INTEGER, PUBLIC ::   ndtadyn = 73       !: Number of dat in one year
50   INTEGER, PUBLIC ::   ndtatot = 73       !: Number of data in the input field
51   INTEGER, PUBLIC ::   nsptint = 1        !: type of spatial interpolation
52
53   CHARACTER(len=45) ::   cfile_grid_T = 'dyna_grid_T.nc'   ! name of the grid_T file
54   CHARACTER(len=45) ::   cfile_grid_U = 'dyna_grid_U.nc'   ! name of the grid_U file
55   CHARACTER(len=45) ::   cfile_grid_V = 'dyna_grid_V.nc'   ! name of the grid_V file
56   CHARACTER(len=45) ::   cfile_grid_W = 'dyna_grid_W.nc'   ! name of the grid_W file
57   
58   REAL(wp) ::   rnspdta    ! number of time step per 2 consecutives data
59   REAL(wp) ::   rnspdta2   ! rnspdta * 0.5
60
61   INTEGER ::   ndyn1, ndyn2    !
62   INTEGER ::   nlecoff = 0     ! switch for the first read
63   INTEGER ::   numfl_t, numfl_u, numfl_v, numfl_w
64
65   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: tdta       ! temperature at two consecutive times
66   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: sdta       ! salinity at two consecutive times
67   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: udta       ! zonal velocity at two consecutive times
68   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: vdta       ! meridional velocity at two consecutive times
69   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: wdta       ! vertical velocity at two consecutive times
70   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: avtdta     ! vertical diffusivity coefficient
71
72   REAL(wp), DIMENSION(jpi,jpj    ,2) :: hmlddta    ! mixed layer depth at two consecutive times
73   REAL(wp), DIMENSION(jpi,jpj    ,2) :: wspddta    ! wind speed at two consecutive times
74   REAL(wp), DIMENSION(jpi,jpj    ,2) :: frlddta    ! sea-ice fraction at two consecutive times
75   REAL(wp), DIMENSION(jpi,jpj    ,2) :: empdta     ! E-P at two consecutive times
76   REAL(wp), DIMENSION(jpi,jpj    ,2) :: qsrdta     ! short wave heat flux at two consecutive times
77   REAL(wp), DIMENSION(jpi,jpj    ,2) :: bblxdta    ! frequency of bbl in the x direction at 2 consecutive times
78   REAL(wp), DIMENSION(jpi,jpj    ,2) :: bblydta    ! frequency of bbl in the y direction at 2 consecutive times
79   LOGICAL :: l_offbbl
80#if defined key_ldfslp
81   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: uslpdta    ! zonal isopycnal slopes
82   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: vslpdta    ! meridional isopycnal slopes
83   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: wslpidta   ! zonal diapycnal slopes
84   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: wslpjdta   ! meridional diapycnal slopes
85#endif
86#if ! defined key_degrad &&  defined key_traldf_c2d && defined key_traldf_eiv
87   REAL(wp), DIMENSION(jpi,jpj    ,2) :: aeiwdta    ! G&M coefficient
88#endif
89#if defined key_degrad
90   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: ahtudta, ahtvdta, ahtwdta   ! Lateral diffusivity
91# if defined key_traldf_eiv
92   REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: aeiudta, aeivdta, aeiwdta   ! G&M coefficient
93# endif
94#endif
95
96   !! * Substitutions
97#  include "domzgr_substitute.h90"
98#  include "vectopt_loop_substitute.h90"
99   !!----------------------------------------------------------------------
100   !! NEMO/OFF 3.3 , NEMO Consortium (2010)
101   !! $Id$
102   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
103   !!----------------------------------------------------------------------
104CONTAINS
105
106   SUBROUTINE dta_dyn( kt )
107      !!----------------------------------------------------------------------
108      !!                  ***  ROUTINE dta_dyn  ***
109      !!
110      !! ** Purpose :   Prepares dynamics and physics fields from an NEMO run
111      !!              for an off-line simulation of passive tracers
112      !!
113      !! ** Method : calculates the position of DATA to read READ DATA
114      !!             (example month changement) computes slopes IF needed
115      !!             interpolates DATA IF needed
116      !!----------------------------------------------------------------------
117      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
118      !!
119      INTEGER  ::   iper, iperm1, iswap, izt   ! local integers
120      REAL(wp) ::   zt, zweigh                 ! local scalars
121      !!----------------------------------------------------------------------
122
123      zt     = ( REAL(kt,wp) + rnspdta2 ) / rnspdta
124      izt    = INT( zt )
125      zweigh = zt - REAL( INT(zt), wp )
126
127      IF( lperdyn ) THEN   ;   iperm1 = MOD( izt, ndtadyn )
128      ELSE                 ;   iperm1 = MOD( izt, ndtatot - 1 ) + 1
129      ENDIF
130
131      iper = iperm1 + 1
132      IF( iperm1 == 0 ) THEN
133          IF( lperdyn ) THEN
134              iperm1 = ndtadyn
135          ELSE
136              IF( lfirdyn ) THEN
137                  IF(lwp) WRITE (numout,*) 'dta_dyn:  dynamic file is not periodic with or without interpolation    &
138                     &                                we take the first value for the previous period iperm1 = 0  '
139              END IF
140          END IF
141      END IF
142
143      iswap  = 0
144
145      ! 1. First call lfirdyn = true
146      ! ----------------------------
147
148      IF( lfirdyn ) THEN
149         ndyn1 = iperm1         ! store the information of the period read
150         ndyn2 = iper
151         
152         IF(lwp) THEN
153            WRITE (numout,*) ' dynamics data read for the period ndyn1 =', ndyn1,   &
154               &             ' and for the period ndyn2 = ', ndyn2
155            WRITE (numout,*) ' time step is : ', kt
156            WRITE (numout,*) ' we have ndtadyn = ', ndtadyn, ' records in the dynamic file for one year'
157         END IF
158         !
159         CALL dynrea( kt, MAX( 1, iperm1) )           ! data read for the iperm1 period
160         
161         IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN      ! Computes slopes (here tsn and avt are used as workspace)
162            tsn (:,:,:,jp_tem) = tdta  (:,:,:,2)
163            tsn (:,:,:,jp_sal) = sdta  (:,:,:,2)
164            avt(:,:,:)         = avtdta(:,:,:,2)
165         
166            CALL eos( tsn, rhd, rhop )   ! Time-filtered in situ density
167            CALL bn2( tsn, rn2 )         ! before Brunt-Vaisala frequency
168            IF( ln_zps )   &
169               &   CALL zps_hde( kt, jpts, tsn, gtsu, gtsv,  &  ! Partial steps: before Horizontal DErivative
170               &                           rhd, gru , grv   )    ! of t, s, rd at the bottom ocean level
171            CALL zdf_mxl( kt )           ! mixed layer depth
172            CALL ldf_slp( kt, rhd, rn2 )
173         
174            uslpdta (:,:,:,2) = uslp (:,:,:)
175            vslpdta (:,:,:,2) = vslp (:,:,:)
176            wslpidta(:,:,:,2) = wslpi(:,:,:)
177            wslpjdta(:,:,:,2) = wslpj(:,:,:)
178         END IF
179         !
180         CALL swap_dyn_data            ! swap from record 2 to 1
181         !
182         iswap = 1        !  indicates swap
183         !
184         CALL dynrea( kt, iper )       ! data read for the iper period
185         !
186         IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN      ! Computes slopes (here tsn and avt are used as workspace)
187            tsn (:,:,:,jp_tem) = tdta  (:,:,:,2)
188            tsn (:,:,:,jp_sal) = sdta  (:,:,:,2)
189            avt(:,:,:)         = avtdta(:,:,:,2)
190            !
191                           CALL eos( tsn, rhd, rhop )                   ! now in situ density
192                           CALL bn2( tsn, rn2 )                         ! now Brunt-Vaisala frequency
193            IF( ln_zps )   CALL zps_hde( kt, jpts, tsn, gtsu, gtsv,  &  ! Partial steps: before Horizontal DErivative
194               &                                   rhd, gru , grv   )   ! of t, s, rd at the bottom ocean level
195                           CALL zdf_mxl( kt )                           ! mixed layer depth
196                           CALL ldf_slp( kt, rhd, rn2 )                 ! slope of iso-neutral surfaces
197            !
198            uslpdta (:,:,:,2) = uslp (:,:,:)
199            vslpdta (:,:,:,2) = vslp (:,:,:)
200            wslpidta(:,:,:,2) = wslpi(:,:,:)
201            wslpjdta(:,:,:,2) = wslpj(:,:,:)
202         END IF
203         !
204         lfirdyn = .FALSE.    ! trace the first call
205      ENDIF
206      !
207      ! And now what we have to do at every time step
208      ! check the validity of the period in memory
209      !
210      IF( iperm1 /= ndyn1 ) THEN 
211         !
212         IF( iperm1 == 0 ) THEN
213            IF(lwp) THEN
214               WRITE (numout,*) ' dynamic file is not periodic with periodic interpolation'
215               WRITE (numout,*) ' we take the last value for the last period '
216               WRITE (numout,*) ' iperm1 = 12,   iper = 13  '
217            ENDIF
218            iperm1 = 12
219            iper   = 13
220         ENDIF
221         !
222         CALL swap_dyn_data         ! We have to prepare a new read of data : swap from record 2 to 1
223         !
224         iswap = 1                  !  indicates swap
225         !
226         CALL dynrea( kt, iper )    ! data read for the iper period
227         !
228         IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN
229            ! Computes slopes. Caution : here tsn and avt are used as workspace
230            tsn(:,:,:,jp_tem) = tdta  (:,:,:,2)
231            tsn(:,:,:,jp_sal) = sdta  (:,:,:,2)
232            avt(:,:,:)        = avtdta(:,:,:,2)
233            !
234                           CALL eos( tsn, rhd, rhop )                   ! now in situ density
235                           CALL bn2( tsn, rn2 )                         ! now Brunt-Vaisala frequency
236            IF( ln_zps )   CALL zps_hde( kt, jpts, tsn, gtsu, gtsv,  &  ! Partial steps: before Horizontal DErivative
237               &                                   rhd, gru , grv   )   ! of t, s, rd at the bottom ocean level
238            CALL zdf_mxl( kt )                                          ! mixed layer depth
239            CALL ldf_slp( kt, rhd, rn2 )                                ! slope of iso-neutral surfaces
240            !
241            uslpdta (:,:,:,2) = uslp (:,:,:)
242            vslpdta (:,:,:,2) = vslp (:,:,:)
243            wslpidta(:,:,:,2) = wslpi(:,:,:)
244            wslpjdta(:,:,:,2) = wslpj(:,:,:)
245         END IF
246         !
247         ndyn1 = ndyn2         ! store the information of the period read
248         ndyn2 = iper
249         !
250         IF(lwp) THEN
251            WRITE (numout,*) ' dynamics data read for the period ndyn1 =', ndyn1,   &
252               &             ' and for the period ndyn2 = ', ndyn2
253            WRITE (numout,*) ' time step is : ', kt
254         END IF
255         !
256      END IF
257      !
258      ! Compute the data at the given time step
259      !----------------------------------------     
260
261      IF( nsptint == 0 ) THEN          ! No space interpolation, data are probably correct
262         !                             ! We have to initialize data if we have changed the period         
263         CALL assign_dyn_data
264      ELSEIF( nsptint == 1 ) THEN      ! linear interpolation
265         CALL linear_interp_dyn_data( zweigh )
266      ELSE                             ! other interpolation
267         WRITE (numout,*) ' this kind of interpolation do not exist at the moment : we stop'
268         STOP 'dtadyn'         
269      END IF
270      !
271      CALL eos( tsn, rhd, rhop )       ! In any case, we need rhop
272      !
273#if ! defined key_degrad && defined key_traldf_c2d
274      !                                ! In case of 2D varying coefficients, we need aeiv and aeiu
275      IF( lk_traldf_eiv )   CALL dta_eiv( kt )      ! eddy induced velocity coefficient
276#endif
277      !
278      IF( .NOT. l_offbbl ) THEN       ! Compute bbl coefficients if needed
279         tsb(:,:,:,:) = tsn(:,:,:,:)
280         CALL bbl( kt, 'TRC')
281      END IF
282      IF(ln_ctl) THEN
283         CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' tn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
284         CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_sal), clinfo1=' sn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
285         CALL prt_ctl(tab3d_1=un               , clinfo1=' un      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
286         CALL prt_ctl(tab3d_1=vn               , clinfo1=' vn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
287         CALL prt_ctl(tab3d_1=wn               , clinfo1=' wn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
288         CALL prt_ctl(tab3d_1=avt              , clinfo1=' kz      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
289         CALL prt_ctl(tab2d_1=fr_i             , clinfo1=' fr_i    - : ', mask1=tmask, ovlap=1 )
290         CALL prt_ctl(tab2d_1=hmld             , clinfo1=' hmld    - : ', mask1=tmask, ovlap=1 )
291         CALL prt_ctl(tab2d_1=emps             , clinfo1=' emps    - : ', mask1=tmask, ovlap=1 )
292         CALL prt_ctl(tab2d_1=wndm             , clinfo1=' wspd    - : ', mask1=tmask, ovlap=1 )
293         CALL prt_ctl(tab2d_1=qsr              , clinfo1=' qsr     - : ', mask1=tmask, ovlap=1 )
294      ENDIF
295      !
296   END SUBROUTINE dta_dyn
297
298
299   SUBROUTINE dynrea( kt, kenr )
300      !!----------------------------------------------------------------------
301      !!                  ***  ROUTINE dynrea  ***
302      !!
303      !! ** Purpose : READ dynamics fiels from OPA9 netcdf output
304      !!
305      !! ** Method : READ the kenr records of DATA and store in udta(...,2), .... 
306      !!----------------------------------------------------------------------
307      INTEGER, INTENT(in) ::   kt, kenr   ! time index
308      !!
309      INTEGER ::  jkenr
310      REAL(wp), DIMENSION(jpi,jpj,jpk) ::  zu, zv, zw, zt, zs, zavt , zhdiv              ! 3D workspace
311      REAL(wp), DIMENSION(jpi,jpj)     ::  zemp, zqsr, zmld, zice, zwspd, ztaux, ztauy   ! 2D workspace
312      REAL(wp), DIMENSION(jpi,jpj)     ::  zbblx, zbbly
313
314#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
315      REAL(wp), DIMENSION(jpi,jpj) :: zaeiw 
316#endif
317#if defined key_degrad
318   REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zahtu, zahtv, zahtw  !  Lateral diffusivity
319# if defined key_traldf_eiv
320   REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zaeiu, zaeiv, zaeiw  ! G&M coefficient
321# endif
322#endif
323      !!----------------------------------------------------------------------
324
325      ! 0. Initialization
326     
327      ! cas d'un fichier non periodique : on utilise deux fois le premier et
328      ! le dernier champ temporel
329
330      jkenr = kenr
331
332      IF(lwp) THEN
333         WRITE(numout,*)
334         WRITE(numout,*) 'Dynrea : read dynamical fields, kenr = ', jkenr
335         WRITE(numout,*) '~~~~~~~'
336#if defined key_degrad
337         WRITE(numout,*) ' Degraded fields'
338#endif
339         WRITE(numout,*)
340      ENDIF
341
342
343      IF( kt == nit000 .AND. nlecoff == 0 ) THEN
344         nlecoff = 1
345         CALL  iom_open ( cfile_grid_T, numfl_t )
346         CALL  iom_open ( cfile_grid_U, numfl_u )
347         CALL  iom_open ( cfile_grid_V, numfl_v )
348         CALL  iom_open ( cfile_grid_W, numfl_w )
349      ENDIF
350
351      ! file grid-T
352      !---------------
353      CALL iom_get( numfl_t, jpdom_data, 'votemper', zt   (:,:,:), jkenr )
354      CALL iom_get( numfl_t, jpdom_data, 'vosaline', zs   (:,:,:), jkenr )
355      CALL iom_get( numfl_t, jpdom_data, 'somixhgt', zmld (:,:  ), jkenr )
356      CALL iom_get( numfl_t, jpdom_data, 'sowaflcd', zemp (:,:  ), jkenr )
357      CALL iom_get( numfl_t, jpdom_data, 'soshfldo', zqsr (:,:  ), jkenr )
358      CALL iom_get( numfl_t, jpdom_data, 'soicecov', zice (:,:  ), jkenr )
359      IF( iom_varid( numfl_t, 'sowindsp', ldstop = .FALSE. ) > 0 ) THEN
360         CALL iom_get( numfl_t, jpdom_data, 'sowindsp', zwspd(:,:), jkenr ) 
361      ELSE
362         CALL iom_get( numfl_u, jpdom_data, 'sozotaux', ztaux(:,:), jkenr )
363         CALL iom_get( numfl_v, jpdom_data, 'sometauy', ztauy(:,:), jkenr )
364         CALL tau2wnd( ztaux, ztauy, zwspd )
365      ENDIF
366      ! files grid-U / grid_V
367      CALL iom_get( numfl_u, jpdom_data, 'vozocrtx', zu   (:,:,:), jkenr )
368      CALL iom_get( numfl_v, jpdom_data, 'vomecrty', zv   (:,:,:), jkenr )
369#if defined key_trabbl
370      IF( .NOT. lk_c1d .AND. nn_bbl_ldf == 1 ) THEN
371         IF( iom_varid( numfl_u, 'sobblcox', ldstop = .FALSE. ) > 0  .AND. &
372         &   iom_varid( numfl_v, 'sobblcoy', ldstop = .FALSE. ) > 0 ) THEN
373             CALL iom_get( numfl_u, jpdom_data, 'sobblcox', zbblx(:,:), jkenr )
374             CALL iom_get( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,:), jkenr )
375             l_offbbl = .TRUE.
376         ENDIF
377      ENDIF
378#endif 
379
380      ! file grid-W
381!!      CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw   (:,:,:), jkenr )
382      ! Computation of vertical velocity using horizontal divergence
383      CALL wzv( zu, zv, zw, zhdiv )
384
385      IF( iom_varid( numfl_w, 'voddmavs', ldstop = .FALSE. ) > 0 ) THEN          ! avs exist: it is used
386         CALL iom_get( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr )
387      ELSE                                                                       ! no avs: use avt
388         CALL iom_get( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr )
389      ENDIF
390
391#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
392      CALL iom_get( numfl_w, jpdom_data, 'soleaeiw', zaeiw (:,: ), jkenr )
393#endif
394
395#if defined key_degrad
396      CALL iom_get( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr )
397      CALL iom_get( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr )
398      CALL iom_get( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr )
399#  if defined key_traldf_eiv
400      CALL iom_get( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr )
401      CALL iom_get( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr )
402      CALL iom_get( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr )
403#  endif
404#endif
405
406      udta(:,:,:,2) = zu(:,:,:) * umask(:,:,:)
407      vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) 
408      wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:)
409
410      tdta(:,:,:,2)   = zt  (:,:,:) * tmask(:,:,:)
411      sdta(:,:,:,2)   = zs  (:,:,:) * tmask(:,:,:)
412      avtdta(:,:,:,2) = zavt(:,:,:) * tmask(:,:,:)
413
414#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
415      aeiwdta(:,:,2)  = zaeiw(:,:) * tmask(:,:,1)
416#endif
417
418#if defined key_degrad
419        ahtudta(:,:,:,2) = zahtu(:,:,:) * umask(:,:,:)
420        ahtvdta(:,:,:,2) = zahtv(:,:,:) * vmask(:,:,:)
421        ahtwdta(:,:,:,2) = zahtw(:,:,:) * tmask(:,:,:)
422#  if defined key_traldf_eiv
423        aeiudta(:,:,:,2) = zaeiu(:,:,:) * umask(:,:,:)
424        aeivdta(:,:,:,2) = zaeiv(:,:,:) * vmask(:,:,:)
425        aeiwdta(:,:,:,2) = zaeiw(:,:,:) * tmask(:,:,:)
426#  endif
427#endif
428
429      ! fluxes
430      !
431      wspddta(:,:,2)  = zwspd(:,:) * tmask(:,:,1)
432      frlddta(:,:,2)  = MIN( 1., zice(:,:) ) * tmask(:,:,1)
433      empdta (:,:,2)  = zemp(:,:) * tmask(:,:,1)
434      qsrdta (:,:,2)  = zqsr(:,:) * tmask(:,:,1)
435      hmlddta(:,:,2)  = zmld(:,:) * tmask(:,:,1)
436
437#if defined key_trabbl
438      IF( l_offbbl ) THEN
439         bblxdta(:,:,2) = MAX( 0., zbblx(:,:) )
440         bblydta(:,:,2) = MAX( 0., zbbly(:,:) )
441         WHERE( bblxdta(:,:,2) > 2. ) bblxdta(:,:,2) = 0.
442         WHERE( bblydta(:,:,2) > 2. ) bblydta(:,:,2) = 0.
443      ENDIF
444#endif
445     
446      IF( kt == nitend ) THEN
447         CALL iom_close ( numfl_t )
448         CALL iom_close ( numfl_u )
449         CALL iom_close ( numfl_v )
450         CALL iom_close ( numfl_w )
451      ENDIF
452      !     
453   END SUBROUTINE dynrea
454
455
456   SUBROUTINE dta_dyn_init
457      !!----------------------------------------------------------------------
458      !!                  ***  ROUTINE dta_dyn_init  ***
459      !!
460      !! ** Purpose :   initializations of parameters for the interpolation
461      !!
462      !! ** Method :
463      !!----------------------------------------------------------------------
464      REAL(wp) ::   znspyr   !: number of time step per year
465      !!
466      NAMELIST/namdyn/ ndtadyn, ndtatot, nsptint, lperdyn,  &
467      &                cfile_grid_T, cfile_grid_U, cfile_grid_V, cfile_grid_W
468      !!----------------------------------------------------------------------
469
470      !  Define the dynamical input parameters
471      ! ======================================
472
473      REWIND( numnam )              ! Read Namelist namdyn : Lateral physics on tracers
474      READ  ( numnam, namdyn )
475
476      IF(lwp) THEN                  ! control print
477         WRITE(numout,*)
478         WRITE(numout,*) 'namdyn : offline dynamical selection'
479         WRITE(numout,*) '~~~~~~~'
480         WRITE(numout,*) '  Namelist namdyn : set parameters for the lecture of the dynamical fields'
481         WRITE(numout,*) 
482         WRITE(numout,*) ' number of elements in the FILE for a year  ndtadyn = ' , ndtadyn
483         WRITE(numout,*) ' total number of elements in the FILE       ndtatot = ' , ndtatot
484         WRITE(numout,*) ' type of interpolation                      nsptint = ' , nsptint
485         WRITE(numout,*) ' loop on the same FILE                      lperdyn = ' , lperdyn
486         WRITE(numout,*) '  '
487         WRITE(numout,*) ' name of grid_T file                   cfile_grid_T = ', TRIM(cfile_grid_T)   
488         WRITE(numout,*) ' name of grid_U file                   cfile_grid_U = ', TRIM(cfile_grid_U) 
489         WRITE(numout,*) ' name of grid_V file                   cfile_grid_V = ', TRIM(cfile_grid_V) 
490         WRITE(numout,*) ' name of grid_W file                   cfile_grid_W = ', TRIM(cfile_grid_W)     
491         WRITE(numout,*) ' '
492      ENDIF
493      !
494      znspyr   = nyear_len(1) * rday / rdt 
495      rnspdta  = znspyr / FLOAT( ndtadyn )
496      rnspdta2 = rnspdta * 0.5 
497      !
498      CALL dta_dyn( nit000 )
499      !
500   END SUBROUTINE dta_dyn_init
501
502
503   SUBROUTINE wzv( pu, pv, pw, phdiv )
504      !!----------------------------------------------------------------------
505      !!                    ***  ROUTINE wzv  ***
506      !!
507      !! ** Purpose :   Compute the now vertical velocity after the array swap
508      !!
509      !! ** Method  : - compute the now divergence given by :
510      !!         * z-coordinate ONLY !!!!
511      !!         hdiv = 1/(e1t*e2t) [ di(e2u  u) + dj(e1v  v) ]
512      !!     - Using the incompressibility hypothesis, the vertical
513      !!      velocity is computed by integrating the horizontal divergence
514      !!      from the bottom to the surface.
515      !!        The boundary conditions are w=0 at the bottom (no flux).
516      !!----------------------------------------------------------------------
517      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) :: pu, pv    !:  horizontal velocities
518      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(  out) :: pw        !:  verticla velocity
519      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) :: phdiv     !:  horizontal divergence
520      !!
521      INTEGER  ::  ji, jj, jk
522      REAL(wp) ::  zu, zu1, zv, zv1, zet
523      !!----------------------------------------------------------------------
524      !
525      ! Computation of vertical velocity using horizontal divergence
526      phdiv(:,:,:) = 0.
527      DO jk = 1, jpkm1
528         DO jj = 2, jpjm1
529            DO ji = fs_2, fs_jpim1   ! vector opt.
530               zu  = pu(ji  ,jj  ,jk) * umask(ji  ,jj  ,jk) * e2u(ji  ,jj  ) * fse3u(ji  ,jj  ,jk)
531               zu1 = pu(ji-1,jj  ,jk) * umask(ji-1,jj  ,jk) * e2u(ji-1,jj  ) * fse3u(ji-1,jj  ,jk)
532               zv  = pv(ji  ,jj  ,jk) * vmask(ji  ,jj  ,jk) * e1v(ji  ,jj  ) * fse3v(ji  ,jj  ,jk)
533               zv1 = pv(ji  ,jj-1,jk) * vmask(ji  ,jj-1,jk) * e1v(ji  ,jj-1) * fse3v(ji  ,jj-1,jk)
534               zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
535               phdiv(ji,jj,jk) = ( zu - zu1 + zv - zv1 ) * zet 
536            END DO
537         END DO
538      END DO
539      CALL lbc_lnk( phdiv, 'T', 1. )      ! Lateral boundary conditions on phdiv
540      !
541      ! computation of vertical velocity from the bottom
542      pw(:,:,jpk) = 0._wp
543      DO jk = jpkm1, 1, -1
544         pw(:,:,jk) = pw(:,:,jk+1) - fse3t(:,:,jk) * phdiv(:,:,jk)
545      END DO
546      !
547   END SUBROUTINE wzv
548
549
550   SUBROUTINE dta_eiv( kt )
551      !!----------------------------------------------------------------------
552      !!                  ***  ROUTINE dta_eiv  ***
553      !!
554      !! ** Purpose :   Compute the eddy induced velocity coefficient from the
555      !!      growth rate of baroclinic instability.
556      !!
557      !! ** Method : Specific to the offline model. Computes the horizontal
558      !!             values from the vertical value
559      !!----------------------------------------------------------------------
560      INTEGER, INTENT( in ) ::   kt     ! ocean time-step inedx
561      !!
562      INTEGER ::   ji, jj           ! dummy loop indices
563      !!----------------------------------------------------------------------
564      !
565      IF( kt == nit000 ) THEN
566         IF(lwp) WRITE(numout,*)
567         IF(lwp) WRITE(numout,*) 'dta_eiv : eddy induced velocity coefficients'
568         IF(lwp) WRITE(numout,*) '~~~~~~~'
569      ENDIF
570      !
571#if defined key_ldfeiv
572      ! Average the diffusive coefficient at u- v- points
573      DO jj = 2, jpjm1
574         DO ji = fs_2, fs_jpim1   ! vector opt.
575            aeiu(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji+1,jj  ) )
576            aeiv(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji  ,jj+1) )
577         END DO
578      END DO
579      CALL lbc_lnk( aeiu, 'U', 1. )   ;   CALL lbc_lnk( aeiv, 'V', 1. )    ! lateral boundary condition
580#endif
581      !
582   END SUBROUTINE dta_eiv
583
584
585   SUBROUTINE tau2wnd( ptaux, ptauy, pwspd )
586      !!---------------------------------------------------------------------
587      !!                    ***  ROUTINE sbc_tau2wnd  ***
588      !!
589      !! ** Purpose : Estimation of wind speed as a function of wind stress
590      !!
591      !! ** Method  : |tau|=rhoa*Cd*|U|^2
592      !!---------------------------------------------------------------------
593      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) ::   ptaux, ptauy   ! wind stress in i-j direction resp.
594      REAL(wp), DIMENSION(jpi,jpj), INTENT(  out) ::   pwspd          ! wind speed
595      !!
596      REAL(wp) ::   zrhoa  = 1.22_wp       ! Air density kg/m3
597      REAL(wp) ::   zcdrag = 1.5e-3_wp     ! drag coefficient
598      REAL(wp) ::   ztx, zty, ztau, zcoef  ! temporary variables
599      INTEGER  ::   ji, jj                 ! dummy indices
600      !!---------------------------------------------------------------------
601      zcoef = 1. / ( zrhoa * zcdrag )
602!CDIR NOVERRCHK
603      DO jj = 2, jpjm1
604!CDIR NOVERRCHK
605         DO ji = fs_2, fs_jpim1   ! vector opt.
606            ztx = ptaux(ji,jj) * umask(ji,jj,1) + ptaux(ji-1,jj  ) * umask(ji-1,jj  ,1)
607            zty = ptauy(ji,jj) * vmask(ji,jj,1) + ptauy(ji  ,jj-1) * vmask(ji  ,jj-1,1)
608            ztau = 0.5 * SQRT( ztx * ztx + zty * zty )
609            pwspd(ji,jj) = SQRT ( ztau * zcoef ) * tmask(ji,jj,1)
610         END DO
611      END DO
612      CALL lbc_lnk( pwspd(:,:), 'T', 1. )
613      !
614   END SUBROUTINE tau2wnd
615
616
617   SUBROUTINE swap_dyn_data
618      !!----------------------------------------------------------------------
619      !!                    ***  ROUTINE swap_dyn_data  ***
620      !!
621      !! ** Purpose :   swap array data
622      !!----------------------------------------------------------------------
623      !
624      ! swap from record 2 to 1
625      tdta   (:,:,:,1) = tdta   (:,:,:,2)
626      sdta   (:,:,:,1) = sdta   (:,:,:,2)
627      avtdta (:,:,:,1) = avtdta (:,:,:,2)
628      udta   (:,:,:,1) = udta   (:,:,:,2)
629      vdta   (:,:,:,1) = vdta   (:,:,:,2)
630      wdta   (:,:,:,1) = wdta   (:,:,:,2)
631#if defined key_ldfslp && ! defined key_c1d
632      uslpdta (:,:,:,1) = uslpdta (:,:,:,2)
633      vslpdta (:,:,:,1) = vslpdta (:,:,:,2)
634      wslpidta(:,:,:,1) = wslpidta(:,:,:,2)
635      wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2)
636#endif
637      hmlddta(:,:,1) = hmlddta(:,:,2) 
638      wspddta(:,:,1) = wspddta(:,:,2) 
639      frlddta(:,:,1) = frlddta(:,:,2) 
640      empdta (:,:,1) = empdta (:,:,2) 
641      qsrdta (:,:,1) = qsrdta (:,:,2) 
642      IF( l_offbbl ) THEN
643         bblxdta(:,:,1) = bblxdta(:,:,2)
644         bblydta(:,:,1) = bblydta(:,:,2) 
645      ENDIF
646
647#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
648      aeiwdta(:,:,1) = aeiwdta(:,:,2)
649#endif
650
651#if defined key_degrad
652      ahtudta(:,:,:,1) = ahtudta(:,:,:,2)
653      ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2)
654      ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2)
655#  if defined key_traldf_eiv
656      aeiudta(:,:,:,1) = aeiudta(:,:,:,2)
657      aeivdta(:,:,:,1) = aeivdta(:,:,:,2)
658      aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2)
659#  endif
660#endif
661      !
662   END SUBROUTINE swap_dyn_data
663
664
665   SUBROUTINE assign_dyn_data
666      !!----------------------------------------------------------------------
667      !!                    ***  ROUTINE assign_dyn_data  ***
668      !!
669      !! ** Purpose :   Assign dynamical data to the data that have been read
670      !!                without time interpolation
671      !!
672      !!----------------------------------------------------------------------
673     
674      tsn(:,:,:,jp_tem) = tdta  (:,:,:,2)
675      tsn(:,:,:,jp_sal) = sdta  (:,:,:,2)
676      avt(:,:,:)        = avtdta(:,:,:,2)
677     
678      un (:,:,:) = udta  (:,:,:,2) 
679      vn (:,:,:) = vdta  (:,:,:,2)
680      wn (:,:,:) = wdta  (:,:,:,2)
681     
682#if defined key_ldfslp && ! defined key_c1d
683      uslp (:,:,:) = uslpdta (:,:,:,2) 
684      vslp (:,:,:) = vslpdta (:,:,:,2) 
685      wslpi(:,:,:) = wslpidta(:,:,:,2) 
686      wslpj(:,:,:) = wslpjdta(:,:,:,2) 
687#endif
688
689      hmld(:,:) = hmlddta(:,:,2) 
690      wndm(:,:) = wspddta(:,:,2) 
691      fr_i(:,:) = frlddta(:,:,2) 
692      emp (:,:) = empdta (:,:,2) 
693      emps(:,:) = emp(:,:) 
694      qsr (:,:) = qsrdta (:,:,2) 
695#if defined key_trabbl
696      IF( l_offbbl ) THEN
697         ahu_bbl(:,:) = bblxdta(:,:,2)
698         ahv_bbl(:,:) = bblydta(:,:,2) 
699      ENDIF
700#endif
701#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
702      aeiw(:,:) = aeiwdta(:,:,2)
703#endif
704     
705#if defined key_degrad
706      ahtu(:,:,:) = ahtudta(:,:,:,2)
707      ahtv(:,:,:) = ahtvdta(:,:,:,2)
708      ahtw(:,:,:) = ahtwdta(:,:,:,2)
709#  if defined key_traldf_eiv
710      aeiu(:,:,:) = aeiudta(:,:,:,2)
711      aeiv(:,:,:) = aeivdta(:,:,:,2)
712      aeiw(:,:,:) = aeiwdta(:,:,:,2)
713#  endif
714#endif
715      !
716   END SUBROUTINE assign_dyn_data
717
718
719   SUBROUTINE linear_interp_dyn_data( pweigh )
720      !!----------------------------------------------------------------------
721      !!               ***  ROUTINE linear_interp_dyn_data  ***
722      !!
723      !! ** Purpose :   linear interpolation of data
724      !!----------------------------------------------------------------------
725      REAL(wp), INTENT(in) ::   pweigh   ! weigh
726      !!
727      REAL(wp) :: zweighm1
728      !!----------------------------------------------------------------------
729
730      zweighm1 = 1. - pweigh
731     
732      tsn(:,:,:,jp_tem) = zweighm1 * tdta  (:,:,:,1) + pweigh * tdta  (:,:,:,2)
733      tsn(:,:,:,jp_sal) = zweighm1 * sdta  (:,:,:,1) + pweigh * sdta  (:,:,:,2)
734      avt(:,:,:)        = zweighm1 * avtdta(:,:,:,1) + pweigh * avtdta(:,:,:,2)
735     
736      un (:,:,:) = zweighm1 * udta  (:,:,:,1) + pweigh * udta  (:,:,:,2) 
737      vn (:,:,:) = zweighm1 * vdta  (:,:,:,1) + pweigh * vdta  (:,:,:,2)
738      wn (:,:,:) = zweighm1 * wdta  (:,:,:,1) + pweigh * wdta  (:,:,:,2)
739     
740#if defined key_ldfslp && ! defined key_c1d
741      uslp (:,:,:) = zweighm1 * uslpdta (:,:,:,1) + pweigh * uslpdta (:,:,:,2) 
742      vslp (:,:,:) = zweighm1 * vslpdta (:,:,:,1) + pweigh * vslpdta (:,:,:,2) 
743      wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + pweigh * wslpidta(:,:,:,2) 
744      wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + pweigh * wslpjdta(:,:,:,2) 
745#endif
746
747      hmld(:,:) = zweighm1 * hmlddta(:,:,1) + pweigh  * hmlddta(:,:,2) 
748      wndm(:,:) = zweighm1 * wspddta(:,:,1) + pweigh  * wspddta(:,:,2) 
749      fr_i(:,:) = zweighm1 * frlddta(:,:,1) + pweigh  * frlddta(:,:,2) 
750      emp (:,:) = zweighm1 * empdta (:,:,1) + pweigh  * empdta (:,:,2) 
751      emps(:,:) = emp(:,:) 
752      qsr (:,:) = zweighm1 * qsrdta (:,:,1) + pweigh  * qsrdta (:,:,2) 
753#if defined key_trabbl
754      IF( l_offbbl ) THEN
755         ahu_bbl(:,:) = zweighm1 * bblxdta(:,:,1) +  pweigh  * bblxdta(:,:,2)
756         ahv_bbl(:,:) = zweighm1 * bblydta(:,:,1) +  pweigh  * bblydta(:,:,2)
757      ENDIF
758#endif
759
760#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
761      aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + pweigh * aeiwdta(:,:,2)
762#endif
763     
764#if defined key_degrad
765      ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + pweigh * ahtudta(:,:,:,2)
766      ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + pweigh * ahtvdta(:,:,:,2)
767      ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + pweigh * ahtwdta(:,:,:,2)
768#  if defined key_traldf_eiv
769      aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + pweigh * aeiudta(:,:,:,2)
770      aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + pweigh * aeivdta(:,:,:,2)
771      aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + pweigh * aeiwdta(:,:,:,2)
772#  endif
773#endif
774      !     
775   END SUBROUTINE linear_interp_dyn_data
776
777   !!======================================================================
778END MODULE dtadyn
Note: See TracBrowser for help on using the repository browser.