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

source: branches/CMIP5_IPSL/NEMO/OFF_SRC/dtadyn.F90 @ 3395

Last change on this file since 3395 was 1873, checked in by cetlod, 14 years ago

Update PISCES diagnostics for AR5

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 40.5 KB
Line 
1MODULE dtadyn
2   !!======================================================================
3   !!                       ***  MODULE  dtadyn  ***
4   !! OFFLINE : interpolation of the physical fields
5   !!=====================================================================
6
7   !!----------------------------------------------------------------------
8   !!   dta_dyn_init : initialization, namelist read, and parameters control
9   !!   dta_dyn      : Interpolation of the fields
10   !!----------------------------------------------------------------------
11   !! * Modules used
12   USE oce             ! ocean dynamics and tracers variables
13   USE dom_oce         ! ocean space and time domain variables
14   USE zdf_oce         ! ocean vertical physics
15   USE in_out_manager  ! I/O manager
16   USE phycst          ! physical constants
17   USE sbc_oce
18   USE ldfslp
19   USE ldfeiv          ! eddy induced velocity coef.      (ldf_eiv routine)
20   USE ldftra_oce      ! ocean tracer   lateral physics
21   USE zdfmxl
22   USE trabbl
23   USE eosbn2
24   USE zdfddm          ! vertical  physics: double diffusion
25   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
26   USE zpshde
27   USE lib_mpp         ! distributed memory computing library
28
29   IMPLICIT NONE
30   PRIVATE
31
32   !! *  Routine accessibility
33   PUBLIC dta_dyn_init   ! called by opa.F90
34   PUBLIC dta_dyn        ! called by step.F90
35
36   LOGICAL , PUBLIC :: &
37      lperdyn = .TRUE. , & ! boolean for periodic fields or not
38      lfirdyn = .TRUE.     ! boolean for the first call or not
39
40   INTEGER , PUBLIC :: &
41      ndtadyn = 73 ,  & ! Number of dat in one year
42      ndtatot = 73 ,  & ! Number of data in the input field
43      nsptint = 1 ,   & ! type of spatial interpolation
44      nficdyn = 2       ! number of dynamical fields
45
46   CHARACTER(len=45)  ::  &
47      cfile_grid_T = 'dyna_grid_T.nc', &   !: name of the grid_T file
48      cfile_grid_U = 'dyna_grid_U.nc', &   !: name of the grid_U file
49      cfile_grid_V = 'dyna_grid_V.nc', &   !: name of the grid_V file
50      cfile_grid_W = 'dyna_grid_W.nc'      !: name of the grid_W file
51   
52   REAL(wp)      ::   &
53      rnspdta  ,       &  !: number of time step per 2 consecutives data
54      rnspdta2            !: rnspdta * 0.5
55
56   INTEGER ::     &
57      ndyn1, ndyn2 , &
58      nlecoff = 0  , & ! switch for the first read
59      numfl_t, numfl_u, &
60      numfl_v, numfl_w
61
62   REAL(wp), DIMENSION(jpi,jpj,jpk,2) ::   &
63      tdta   ,   & ! temperature at two consecutive times
64      sdta   ,   & ! salinity at two consecutive times
65      udta   ,   & ! zonal velocity at two consecutive times
66      vdta   ,   & ! meridional velocity at two consecutive times
67      wdta   ,   & ! vertical velocity at two consecutive times
68#if defined key_trc_diatrd
69      hdivdta,   & ! horizontal divergence
70#endif
71      avtdta       ! vertical diffusivity coefficient
72
73   REAL(wp), DIMENSION(jpi,jpj,2) ::       &
74      hmlddta,   & ! mixed layer depth at two consecutive times
75      wspddta,   & ! wind speed at two consecutive times
76      frlddta,   & ! sea-ice fraction at two consecutive times
77      empdta ,   & ! E-P at two consecutive times
78      qsrdta       ! short wave heat flux at two consecutive times
79
80#if defined key_ldfslp
81   REAL(wp), DIMENSION(jpi,jpj,jpk,2) ::   &
82      uslpdta ,  & ! zonal isopycnal slopes
83      vslpdta ,  & ! meridional isopycnal slopes
84      wslpidta , & ! zonal diapycnal slopes
85      wslpjdta     ! meridional diapycnal slopes
86#endif
87
88#if ! defined key_off_degrad &&  defined key_traldf_c2d
89   REAL(wp), DIMENSION(jpi,jpj,2) ::   &
90      ahtwdta      ! Lateral diffusivity
91# if defined key_trcldf_eiv 
92   REAL(wp), DIMENSION(jpi,jpj,2) ::   &
93      aeiwdta      ! G&M coefficient
94# endif
95#endif
96
97#if defined key_off_degrad
98   REAL(wp), DIMENSION(jpi,jpj,jpk,2) ::   &
99      ahtudta, ahtvdta, ahtwdta  !  Lateral diffusivity
100# if defined key_trcldf_eiv
101   REAL(wp), DIMENSION(jpi,jpj,jpk,2) ::   &
102      aeiudta, aeivdta, aeiwdta  ! G&M coefficient
103# endif
104
105#endif
106
107#if defined key_trcbbl_dif   ||   defined key_trcbbl_adv
108   REAL(wp), DIMENSION(jpi,jpj,2) ::       &
109      bblxdta ,  & ! frequency of bbl in the x direction at 2 consecutive times
110      bblydta      ! frequency of bbl in the y direction at 2 consecutive times
111#endif
112
113   !! * Substitutions
114#  include "domzgr_substitute.h90"
115#  include "vectopt_loop_substitute.h90"
116   !!----------------------------------------------------------------------
117   !!   OPA 9.0 , LOCEAN-IPSL  (2005)
118   !!   $Id$
119   !!   This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
120   !!----------------------------------------------------------------------
121
122CONTAINS
123
124   SUBROUTINE dta_dyn( kt )
125      !!----------------------------------------------------------------------
126      !!                  ***  ROUTINE dta_dyn  ***
127      !!
128      !! ** Purpose : Prepares dynamics and physics fields from an
129      !!              OPA9 simulation  for an off-line simulation
130      !!               for passive tracer
131      !!
132      !! ** Method : calculates the position of DATA to read READ DATA
133      !!             (example month changement) computes slopes IF needed
134      !!             interpolates DATA IF needed
135      !!
136      !! ** History :
137      !!   ! original  : 92-01 (M. Imbard: sub domain)
138      !!   ! addition  : 98-04 (L.Bopp MA Foujols: slopes for isopyc.)
139      !!   ! addition  : 98-05 (L. Bopp read output of coupled run)
140      !!   ! addition  : 05-03 (O. Aumont and A. El Moussaoui) F90
141      !!   ! addition  : 05-12 (C. Ethe) Adapted for DEGINT
142      !!----------------------------------------------------------------------
143      !! * Arguments
144      INTEGER, INTENT( in ) ::   kt       ! ocean time-step index
145
146      !! * Local declarations
147      INTEGER ::   iper, iperm1, iswap, izt   
148
149      REAL(wp) :: zt, zweigh
150
151      ! 0. Initialization
152      ! -----------------
153
154      IF( lfirdyn ) THEN
155         ! first time step MUST BE nit000
156         IF( kt /= nit000 ) THEN
157            IF (lwp) THEN
158               WRITE (numout,*) ' kt MUST BE EQUAL to nit000. kt = ',kt ,' nit000 = ',nit000 
159              STOP 'dtadyn'
160            ENDIF
161          ENDIF 
162          ! Initialize the parameters of the interpolation
163          CALL dta_dyn_init
164      ENDIF
165
166      zt       = ( FLOAT (kt) + rnspdta2 ) / rnspdta
167      izt      = INT( zt )
168      zweigh   = zt - FLOAT( INT(zt) )
169
170      IF( lperdyn ) THEN
171         iperm1 = MOD( izt, ndtadyn )
172      ELSE
173         iperm1 = MOD( izt, ndtatot - 1 ) + 1
174      ENDIF
175
176      iper = iperm1 + 1
177      IF( iperm1 == 0 ) THEN
178          IF( lperdyn ) THEN
179              iperm1 = ndtadyn
180          ELSE
181              IF( lfirdyn ) THEN
182                  IF (lwp) WRITE (numout,*) & 
183                      &   ' dynamic file is not periodic with or without interpolation  &
184                      &   we take the first value for the previous period iperm1 = 0  '
185              END IF
186          END IF
187      END IF
188
189      iswap  = 0
190
191      ! 1. First call lfirdyn = true
192      ! ----------------------------
193
194      IF( lfirdyn ) THEN
195         ! store the information of the period read
196         ndyn1 = iperm1
197         ndyn2 = iper
198         
199         IF (lwp) THEN
200            WRITE (numout,*) ' dynamics data read for the period ndyn1 =',ndyn1, &
201               &             ' and for the period ndyn2 = ',ndyn2
202            WRITE (numout,*) ' time step is : ', kt
203            WRITE (numout,*) ' we have ndtadyn = ',ndtadyn,' records in the dynamic file for one year'
204         END IF
205         !
206         IF( iperm1 /= 0 ) THEN         ! data read for the iperm1 period
207            CALL dynrea( kt, iperm1 ) 
208         ELSE
209            CALL dynrea( kt, 1 )
210         ENDIF
211         
212#if defined key_ldfslp
213         ! Computes slopes
214         ! Caution : here tn, sn and avt are used as workspace
215         tn (:,:,:) = tdta  (:,:,:,2)
216         sn (:,:,:) = sdta  (:,:,:,2)
217         avt(:,:,:) = avtdta(:,:,:,2)
218         
219         CALL eos( tn, sn, rhd, rhop )   ! Time-filtered in situ density
220         CALL bn2( tn, sn, rn2 )         ! before Brunt-Vaisala frequency
221         IF( ln_zps )   &
222            &   CALL zps_hde( kt, tn , sn , rhd,  &  ! Partial steps: before Horizontal DErivative
223            &                 gtu, gsu, gru,  &  ! of t, s, rd at the bottom ocean level
224            &                 gtv, gsv, grv )
225         CALL zdf_mxl( kt )              ! mixed layer depth
226         CALL ldf_slp( kt, rhd, rn2 )
227         
228         uslpdta (:,:,:,2) = uslp (:,:,:)
229         vslpdta (:,:,:,2) = vslp (:,:,:)
230         wslpidta(:,:,:,2) = wslpi(:,:,:)
231         wslpjdta(:,:,:,2) = wslpj(:,:,:)
232#endif
233         
234         ! swap from record 2 to 1
235         CALL swap_dyn_data
236         
237         iswap = 1        !  indicates swap
238         
239         CALL dynrea( kt, iper )    ! data read for the iper period
240         
241#if defined key_ldfslp
242         ! Computes slopes
243         ! Caution : here tn, sn and avt are used as workspace
244         tn (:,:,:) = tdta  (:,:,:,2)
245         sn (:,:,:) = sdta  (:,:,:,2)
246         avt(:,:,:) = avtdta(:,:,:,2)
247         
248         CALL eos( tn, sn, rhd, rhop )   ! Time-filtered in situ density
249         CALL bn2( tn, sn, rn2 )         ! before Brunt-Vaisala frequency
250         IF( ln_zps )   &
251            &   CALL zps_hde( kt, tn , sn , rhd,  &  ! Partial steps: before Horizontal DErivative
252            &                 gtu, gsu, gru,  &  ! of t, s, rd at the bottom ocean level
253            &                 gtv, gsv, grv )
254         CALL zdf_mxl( kt )              ! mixed layer depth
255         CALL ldf_slp( kt, rhd, rn2 )
256         
257         uslpdta (:,:,:,2) = uslp (:,:,:)
258         vslpdta (:,:,:,2) = vslp (:,:,:)
259         wslpidta(:,:,:,2) = wslpi(:,:,:)
260         wslpjdta(:,:,:,2) = wslpj(:,:,:)
261#endif
262         !
263         lfirdyn=.FALSE.    ! trace the first call
264      ENDIF
265      !
266      ! And now what we have to do at every time step
267      ! check the validity of the period in memory
268      !
269      IF( iperm1 /= ndyn1 ) THEN
270
271         IF( iperm1 == 0. ) THEN
272            IF (lwp) THEN
273               WRITE (numout,*) ' dynamic file is not periodic with periodic interpolation'
274               WRITE (numout,*) ' we take the last value for the last period '
275               WRITE (numout,*) ' iperm1 = 12,   iper = 13  '
276            ENDIF
277            iperm1 = 12
278            iper   = 13
279         ENDIF
280         !
281         ! We have to prepare a new read of data : swap from record 2 to 1
282         !
283         CALL swap_dyn_data
284
285         iswap = 1        !  indicates swap
286         
287         CALL dynrea( kt, iper )    ! data read for the iper period
288
289#if defined key_ldfslp
290         ! Computes slopes
291         ! Caution : here tn, sn and avt are used as workspace
292         tn (:,:,:) = tdta  (:,:,:,2)
293         sn (:,:,:) = sdta  (:,:,:,2)
294         avt(:,:,:) = avtdta(:,:,:,2)
295         
296         CALL eos( tn, sn, rhd, rhop )   ! Time-filtered in situ density
297         CALL bn2( tn, sn, rn2 )         ! before Brunt-Vaisala frequency
298         IF( ln_zps )   &
299            &   CALL zps_hde( kt, tn , sn , rhd,  &  ! Partial steps: before Horizontal DErivative
300            &                 gtu, gsu, gru,  &  ! of t, s, rd at the bottom ocean level
301            &                 gtv, gsv, grv )
302         CALL zdf_mxl( kt )              ! mixed layer depth
303         CALL ldf_slp( kt, rhd, rn2 )
304         
305         uslpdta (:,:,:,2) = uslp (:,:,:)
306         vslpdta (:,:,:,2) = vslp (:,:,:)
307         wslpidta(:,:,:,2) = wslpi(:,:,:)
308         wslpjdta(:,:,:,2) = wslpj(:,:,:)
309#endif
310       
311         ! store the information of the period read
312         ndyn1 = ndyn2
313         ndyn2 = iper
314
315         IF (lwp) THEN
316            WRITE (numout,*) ' dynamics data read for the period ndyn1 =',ndyn1, &
317               &             ' and for the period ndyn2 = ',ndyn2
318            WRITE (numout,*) ' time step is : ', kt
319         END IF
320         !
321      END IF
322      !
323      ! Compute the data at the given time step
324      !----------------------------------------     
325
326      IF( nsptint == 0 ) THEN   
327         ! No spatial interpolation, data are probably correct
328         ! We have to initialize data if we have changed the period         
329         CALL assign_dyn_data         
330      ELSE IF( nsptint == 1 ) THEN
331         ! linear interpolation
332         CALL linear_interp_dyn_data( zweigh )
333      ELSE 
334         ! other interpolation
335         WRITE (numout,*) ' this kind of interpolation do not exist at the moment : we stop'
336         STOP 'dtadyn'         
337      END IF
338     
339      ! In any case, we need rhop
340      CALL eos( tn, sn, rhd, rhop ) 
341     
342#if ! defined key_off_degrad && defined key_traldf_c2d
343      ! In case of 2D varying coefficients, we need aeiv and aeiu
344      IF( lk_traldf_eiv )   CALL ldf_eiv( kt )      ! eddy induced velocity coefficient
345#endif
346   
347   END SUBROUTINE dta_dyn
348
349   SUBROUTINE dynrea( kt, kenr )
350      !!----------------------------------------------------------------------
351      !!                  ***  ROUTINE dynrea  ***
352      !!
353      !! ** Purpose : READ dynamics fiels from OPA9 netcdf output
354      !!
355      !! ** Method : READ the kenr records of DATA and store in
356      !!             in udta(...,2), .... 
357      !!
358      !! ** History : additions : M. Levy et M. Benjelloul jan 2001
359      !!              (netcdf FORMAT)
360      !!              05-03 (O. Aumont and A. El Moussaoui) F90
361      !!              06-07 : (C. Ethe) use of iom module
362      !!----------------------------------------------------------------------
363      !! * Modules used
364      USE iom
365
366      !! * Arguments
367      INTEGER, INTENT( in ) ::   kt, kenr       ! time index
368      !! * Local declarations
369      INTEGER ::  jkenr
370
371      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   &
372        zu, zv, zw, zt, zs, zavt ,   &     ! 3-D dynamical fields
373        zhdiv                              ! horizontal divergence
374
375      REAL(wp), DIMENSION(jpi,jpj) :: &
376         zemp, zqsr, zmld, zice, zwspd, &
377         ztaux, ztauy
378#if defined key_trcbbl_dif   ||   defined key_trcbbl_adv
379      REAL(wp), DIMENSION(jpi,jpj) :: zbblx, zbbly
380#endif
381
382#if ! defined key_off_degrad && defined key_traldf_c2d
383      REAL(wp), DIMENSION(jpi,jpj) :: zahtw 
384#   if defined key_trcldf_eiv
385      REAL(wp), DIMENSION(jpi,jpj) :: zaeiw 
386#  endif
387#endif
388
389#if defined key_off_degrad
390   REAL(wp), DIMENSION(jpi,jpj,jpk) ::   &
391      zahtu, zahtv, zahtw  !  Lateral diffusivity
392# if defined key_trcldf_eiv
393   REAL(wp), DIMENSION(jpi,jpj,jpk) ::   &
394      zaeiu, zaeiv, zaeiw  ! G&M coefficient
395# endif
396#endif
397
398      !---------------------------------------------------------------
399      ! 0. Initialization
400     
401      ! cas d'un fichier non periodique : on utilise deux fois le premier et
402      ! le dernier champ temporel
403
404      jkenr = kenr
405
406      IF(lwp) THEN
407         WRITE(numout,*)
408         WRITE(numout,*) 'Dynrea : reading dynamical fields, kenr = ', jkenr
409         WRITE(numout,*) ' ~~~~~~~'
410#if defined key_off_degrad
411         WRITE(numout,*) ' Degraded fields'
412#endif
413         WRITE(numout,*)
414      ENDIF
415
416
417      IF( kt == nit000 .AND. nlecoff == 0 ) THEN
418         nlecoff = 1
419         CALL  iom_open ( cfile_grid_T, numfl_t )
420         CALL  iom_open ( cfile_grid_U, numfl_u )
421         CALL  iom_open ( cfile_grid_V, numfl_v )
422         CALL  iom_open ( cfile_grid_W, numfl_w )
423      ENDIF
424
425      ! file grid-T
426      !---------------
427      CALL iom_get( numfl_t, jpdom_data, 'votemper', zt   (:,:,:), jkenr )
428      CALL iom_get( numfl_t, jpdom_data, 'vosaline', zs   (:,:,:), jkenr )
429      CALL iom_get( numfl_t, jpdom_data, 'somixhgt', zmld (:,:  ), jkenr )
430      CALL iom_get( numfl_t, jpdom_data, 'sowaflcd', zemp (:,:  ), jkenr )
431      CALL iom_get( numfl_t, jpdom_data, 'soshfldo', zqsr (:,:  ), jkenr )
432      CALL iom_get( numfl_t, jpdom_data, 'soicecov', zice (:,:  ), jkenr )
433      IF( iom_varid( numfl_t, 'sowindsp', ldstop = .FALSE. ) > 0 ) THEN
434         CALL iom_get( numfl_t, jpdom_data, 'sowindsp', zwspd(:,:), jkenr ) 
435      ELSE
436         CALL iom_get( numfl_u, jpdom_data, 'sozotaux', ztaux(:,:), jkenr )
437         CALL iom_get( numfl_v, jpdom_data, 'sometauy', ztauy(:,:), jkenr )
438         CALL tau2wnd( ztaux, ztauy, zwspd )
439      ENDIF
440      ! files grid-U / grid_V
441      CALL iom_get( numfl_u, jpdom_data, 'vozocrtx', zu   (:,:,:), jkenr )
442      CALL iom_get( numfl_v, jpdom_data, 'vomecrty', zv   (:,:,:), jkenr )
443
444#if defined key_trcbbl_dif || defined key_trcbbl_adv
445      IF( iom_varid( numfl_u, 'sobblcox', ldstop = .FALSE. ) > 0  .AND. &
446      &   iom_varid( numfl_v, 'sobblcoy', ldstop = .FALSE. ) > 0 ) THEN
447         CALL iom_get( numfl_u, jpdom_data, 'sobblcox', zbblx(:,:), jkenr )
448         CALL iom_get( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,:), jkenr )
449      ELSE
450         CALL bbl_sign( zt, zs, zbblx, zbbly )   
451      ENDIF
452#endif
453
454      ! file grid-W
455!!      CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw   (:,:,:), jkenr )
456      ! Computation of vertical velocity using horizontal divergence
457      CALL wzv( zu, zv, zw, zhdiv )
458
459# if defined key_zdfddm
460      CALL iom_get( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr )
461#else
462      CALL iom_get( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr )
463#endif 
464
465#if ! defined key_off_degrad && defined key_traldf_c2d
466      CALL iom_get( numfl_w, jpdom_data, 'soleahtw', zahtw (:,: ), jkenr )
467#  if   defined key_trcldf_eiv 
468      CALL iom_get( numfl_w, jpdom_data, 'soleaeiw', zaeiw (:,: ), jkenr )
469#  endif
470#endif
471
472#if defined key_off_degrad
473      CALL iom_get( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr )
474      CALL iom_get( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr )
475      CALL iom_get( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr )
476#  if defined key_trcldf_eiv
477      CALL iom_get( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr )
478      CALL iom_get( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr )
479      CALL iom_get( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr )
480#  endif
481#endif
482
483      udta(:,:,:,2) = zu(:,:,:) * umask(:,:,:)
484      vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) 
485      wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:)
486
487#if defined key_trc_diatrd
488      hdivdta(:,:,:,2) = zhdiv(:,:,:) * tmask(:,:,:)
489#endif
490
491      tdta(:,:,:,2)   = zt  (:,:,:) * tmask(:,:,:)
492      sdta(:,:,:,2)   = zs  (:,:,:) * tmask(:,:,:)
493      avtdta(:,:,:,2) = zavt(:,:,:) * tmask(:,:,:)
494
495#if ! defined key_off_degrad && defined key_traldf_c2d
496      ahtwdta(:,:,2)  = zahtw(:,:) * tmask(:,:,1)
497#if defined key_trcldf_eiv
498      aeiwdta(:,:,2)  = zaeiw(:,:) * tmask(:,:,1)
499#endif
500#endif
501
502#if defined key_off_degrad
503        ahtudta(:,:,:,2) = zahtu(:,:,:) * umask(:,:,:)
504        ahtvdta(:,:,:,2) = zahtv(:,:,:) * vmask(:,:,:)
505        ahtwdta(:,:,:,2) = zahtw(:,:,:) * tmask(:,:,:)
506#  if defined key_trcldf_eiv
507        aeiudta(:,:,:,2) = zaeiu(:,:,:) * umask(:,:,:)
508        aeivdta(:,:,:,2) = zaeiv(:,:,:) * vmask(:,:,:)
509        aeiwdta(:,:,:,2) = zaeiw(:,:,:) * tmask(:,:,:)
510#  endif
511#endif
512
513      ! fluxes
514      !
515      wspddta(:,:,2)  = zwspd(:,:) * tmask(:,:,1)
516      frlddta(:,:,2)  = MIN( 1., zice(:,:) ) * tmask(:,:,1)
517      empdta (:,:,2)  = zemp(:,:) * tmask(:,:,1)
518      qsrdta (:,:,2)  = zqsr(:,:) * tmask(:,:,1)
519      hmlddta(:,:,2)  = zmld(:,:) * tmask(:,:,1)
520     
521#if defined key_trcbbl_dif   ||   defined key_trcbbl_adv
522      bblxdta(:,:,2) = MAX( 0., zbblx(:,:) )
523      bblydta(:,:,2) = MAX( 0., zbbly(:,:) )
524
525      WHERE( bblxdta(:,:,2) > 2. ) bblxdta(:,:,2) = 0.
526      WHERE( bblydta(:,:,2) > 2. ) bblydta(:,:,2) = 0.
527#endif
528
529      IF( kt == nitend ) THEN
530         CALL iom_close ( numfl_t )
531         CALL iom_close ( numfl_u )
532         CALL iom_close ( numfl_v )
533         CALL iom_close ( numfl_w )
534      ENDIF
535     
536   END SUBROUTINE dynrea
537
538   SUBROUTINE dta_dyn_init
539      !!----------------------------------------------------------------------
540      !!                  ***  ROUTINE dta_dyn_init  ***
541      !!
542      !! ** Purpose :   initializations of parameters for the interpolation
543      !!
544      !! ** Method :
545      !!
546      !! History :
547      !!    ! original  : 92-01 (M. Imbard: sub domain)
548      !!    ! 98-04 (L.Bopp MA Foujols: slopes for isopyc.)
549      !!    ! 98-05 (L. Bopp read output of coupled run)
550      !!    ! 05-03 (O. Aumont and A. El Moussaoui) F90
551      !!----------------------------------------------------------------------
552      !! * Modules used
553
554      !! * Local declarations
555
556      REAL(wp) ::   znspyr   !: number of time step per year
557
558      NAMELIST/namdyn/ ndtadyn, ndtatot, nsptint, nficdyn, lperdyn,  &
559      &                cfile_grid_T, cfile_grid_U, cfile_grid_V, cfile_grid_W
560      !!----------------------------------------------------------------------
561
562      !  Define the dynamical input parameters
563      ! ======================================
564
565      ! Read Namelist namdyn : Lateral physics on tracers
566      REWIND( numnam )
567      READ  ( numnam, namdyn )
568
569      IF(lwp) THEN
570         WRITE(numout,*)
571         WRITE(numout,*) 'namdyn : offline dynamical selection'
572         WRITE(numout,*) '~~~~~~~'
573         WRITE(numout,*) '  Namelist namdyn : set parameters for the lecture of the dynamical fields'
574         WRITE(numout,*) 
575         WRITE(numout,*) ' number of elements in the FILE for a year  ndtadyn = ' , ndtadyn
576         WRITE(numout,*) ' total number of elements in the FILE       ndtatot = ' , ndtatot
577         WRITE(numout,*) ' type of interpolation                      nsptint = ' , nsptint
578         WRITE(numout,*) ' number of dynamics FILE                    nficdyn = ' , nficdyn
579         WRITE(numout,*) ' loop on the same FILE                      lperdyn = ' , lperdyn
580         WRITE(numout,*) '  '
581         WRITE(numout,*) ' name of grid_T file                   cfile_grid_T = ', TRIM(cfile_grid_T)   
582         WRITE(numout,*) ' name of grid_U file                   cfile_grid_U = ', TRIM(cfile_grid_U) 
583         WRITE(numout,*) ' name of grid_V file                   cfile_grid_V = ', TRIM(cfile_grid_V) 
584         WRITE(numout,*) ' name of grid_W file                   cfile_grid_W = ', TRIM(cfile_grid_W)     
585         WRITE(numout,*) ' '
586      ENDIF
587
588      znspyr   = nyear_len(1) * rday / rdt 
589      rnspdta  = znspyr / FLOAT( ndtadyn )
590      rnspdta2 = rnspdta * 0.5 
591
592   END SUBROUTINE dta_dyn_init
593
594   SUBROUTINE wzv( pu, pv, pw, phdiv )
595      !!----------------------------------------------------------------------
596      !!                    ***  ROUTINE wzv  ***
597      !!
598      !! ** Purpose :   Compute the now vertical velocity after the array swap
599      !!
600      !! ** Method  :
601      !! ** Method  : - Divergence:
602      !!      - compute the now divergence given by :
603      !!         * z-coordinate
604      !!         hdiv = 1/(e1t*e2t) [ di(e2u  u) + dj(e1v  v) ]
605      !!     - Using the incompressibility hypothesis, the vertical
606      !!      velocity is computed by integrating the horizontal divergence
607      !!      from the bottom to the surface.
608      !!        The boundary conditions are w=0 at the bottom (no flux) and,
609      !!      in regid-lid case, w=0 at the sea surface.
610      !!
611      !!
612      !! History :
613      !!   9.0  !  02-07  (G. Madec)  Vector optimization
614      !!----------------------------------------------------------------------
615      !! * Arguments
616      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in  )   :: pu, pv    !:  horizontal velocities
617      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( out )   :: pw        !:  verticla velocity
618      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( inout ) :: phdiv     !:  horizontal divergence
619
620      !! * Local declarations
621      INTEGER  ::  ji, jj, jk
622      REAL(wp) ::  zu, zu1, zv, zv1, zet
623
624
625      ! Computation of vertical velocity using horizontal divergence
626      phdiv(:,:,:) = 0.
627      DO jk = 1, jpkm1
628         DO jj = 2, jpjm1
629            DO ji = fs_2, fs_jpim1   ! vector opt.
630#if defined key_zco
631               zu  = pu(ji  ,jj  ,jk) * umask(ji  ,jj  ,jk) * e2u(ji  ,jj  )
632               zu1 = pu(ji-1,jj  ,jk) * umask(ji-1,jj  ,jk) * e2u(ji-1,jj  )
633               zv  = pv(ji  ,jj  ,jk) * vmask(ji  ,jj  ,jk) * e1v(ji  ,jj  )
634               zv1 = pv(ji  ,jj-1,jk) * vmask(ji  ,jj-1,jk) * e1v(ji  ,jj-1)
635               zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) )
636#else
637               zu  = pu(ji  ,jj  ,jk) * umask(ji  ,jj  ,jk) * e2u(ji  ,jj  ) * fse3u(ji  ,jj  ,jk)
638               zu1 = pu(ji-1,jj  ,jk) * umask(ji-1,jj  ,jk) * e2u(ji-1,jj  ) * fse3u(ji-1,jj  ,jk)
639               zv  = pv(ji  ,jj  ,jk) * vmask(ji  ,jj  ,jk) * e1v(ji  ,jj  ) * fse3v(ji  ,jj  ,jk)
640               zv1 = pv(ji  ,jj-1,jk) * vmask(ji  ,jj-1,jk) * e1v(ji  ,jj-1) * fse3v(ji  ,jj-1,jk)
641               zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
642#endif
643               phdiv(ji,jj,jk) = ( zu - zu1 + zv - zv1 ) * zet 
644            END DO
645         END DO
646      ENDDO
647
648      ! Lateral boundary conditions on phdiv
649      CALL lbc_lnk( phdiv, 'T', 1. )
650
651
652      ! computation of vertical velocity from the bottom
653      pw(:,:,jpk) = 0.
654      DO jk = jpkm1, 1, -1
655         pw(:,:,jk) = pw(:,:,jk+1) - fse3t(:,:,jk) * phdiv(:,:,jk)
656      END DO
657
658   END SUBROUTINE wzv
659
660   SUBROUTINE tau2wnd( ptaux, ptauy, pwspd )
661      !!---------------------------------------------------------------------
662      !!                    ***  ROUTINE sbc_tau2wnd  ***
663      !!
664      !! ** Purpose : Estimation of wind speed as a function of wind stress
665      !!
666      !! ** Method  : |tau|=rhoa*Cd*|U|^2
667      !!---------------------------------------------------------------------
668      !! * Arguments
669      REAL(wp), DIMENSION(jpi,jpj), INTENT( in  ) ::  &
670         ptaux, ptauy                              !: wind stress in i-j direction resp.
671      REAL(wp), DIMENSION(jpi,jpj), INTENT( out ) ::  &
672         pwspd                                     !: wind speed
673      REAL(wp) ::   zrhoa  = 1.22         ! Air density kg/m3
674      REAL(wp) ::   zcdrag = 1.5e-3       ! drag coefficient
675      REAL(wp) ::   ztx, zty, ztau, zcoef ! temporary variables
676      INTEGER  ::   ji, jj                ! dummy indices
677      !!---------------------------------------------------------------------
678      zcoef = 1. / ( zrhoa * zcdrag )
679!CDIR NOVERRCHK
680      DO jj = 2, jpjm1
681!CDIR NOVERRCHK
682         DO ji = fs_2, fs_jpim1   ! vector opt.
683            ztx = ptaux(ji,jj) * umask(ji,jj,1) + ptaux(ji-1,jj  ) * umask(ji-1,jj  ,1)
684            zty = ptauy(ji,jj) * vmask(ji,jj,1) + ptauy(ji  ,jj-1) * vmask(ji  ,jj-1,1)
685            ztau = 0.5 * SQRT( ztx * ztx + zty * zty )
686            pwspd(ji,jj) = SQRT ( ztau * zcoef ) * tmask(ji,jj,1)
687         END DO
688      END DO
689      CALL lbc_lnk( pwspd(:,:), 'T', 1. )
690
691   END SUBROUTINE tau2wnd
692
693#if defined key_trcbbl_dif   ||   defined key_trcbbl_adv
694
695   SUBROUTINE bbl_sign( ptn, psn, pbblx, pbbly )
696      !!----------------------------------------------------------------------
697      !!                    ***  ROUTINE bbl_sign  ***
698      !!
699      !! ** Purpose :   Compute the sign of local gradient of density multiplied by the slope
700      !!                along the bottom slope gradient : grad( rho) * grad(h)
701      !!                Need to compute the diffusive bottom boundary layer
702      !!
703      !! ** Method  :   When the product grad( rho) * grad(h) < 0 (where grad
704      !!      is an along bottom slope gradient) an additional lateral diffu-
705      !!      sive trend along the bottom slope is added to the general tracer
706      !!      trend, otherwise nothing is done. See trcbbl.F90
707      !!
708      !!
709      !! History :
710      !!   9.0  !  02-07  (G. Madec)  Vector optimization
711      !!----------------------------------------------------------------------
712      !! * Arguments
713      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in  ) ::  &
714         ptn             ,  &                           !: temperature
715         psn                                            !: salinity
716      REAL(wp), DIMENSION(jpi,jpj), INTENT( out ) ::  &
717         pbblx , pbbly                                  !: sign of bbl in i-j direction resp.
718     
719      !! * Local declarations
720      INTEGER  ::   ji, jj                   ! dummy loop indices
721      INTEGER  ::   ik
722      REAL(wp) ::   &
723         ztx, zsx, zhx, zalbetx, zgdrhox,     &  ! temporary scalars
724         zty, zsy, zhy, zalbety, zgdrhoy 
725      REAL(wp), DIMENSION(jpi,jpj) ::    &
726        ztnb, zsnb, zdep
727      REAL(wp) ::    fsalbt, pft, pfs, pfh   ! statement function
728      !!----------------------------------------------------------------------
729      ! ratio alpha/beta
730      ! ================
731      !  fsalbt: ratio of thermal over saline expension coefficients
732      !       pft :  potential temperature in degrees celcius
733      !       pfs :  salinity anomaly (s-35) in psu
734      !       pfh :  depth in meters
735
736      fsalbt( pft, pfs, pfh ) =                                              &
737         ( ( ( -0.255019e-07 * pft + 0.298357e-05 ) * pft                    &
738                                   - 0.203814e-03 ) * pft                    &
739                                   + 0.170907e-01 ) * pft                    &
740                                   + 0.665157e-01                            &
741         +(-0.678662e-05 * pfs - 0.846960e-04 * pft + 0.378110e-02 ) * pfs   &
742         +  ( ( - 0.302285e-13 * pfh                                         &
743                - 0.251520e-11 * pfs                                         &
744                + 0.512857e-12 * pft * pft          ) * pfh                  &
745                                     - 0.164759e-06   * pfs                  &
746             +(   0.791325e-08 * pft - 0.933746e-06 ) * pft                  &
747                                     + 0.380374e-04 ) * pfh
748
749      ! 0. 2D fields of bottom temperature and salinity, and bottom slope
750      ! -----------------------------------------------------------------
751      ! mbathy= number of w-level, minimum value=1 (cf domrea.F90)
752#  if defined key_vectopt_loop
753      jj = 1
754      DO ji = 1, jpij   ! vector opt. (forced unrolling)
755#  else
756      DO jj = 1, jpj
757         DO ji = 1, jpi
758#  endif
759            ik          =  MAX( mbathy(ji,jj) - 1, 1 )    ! vertical index of the bottom ocean T-level
760            ztnb(ji,jj) = ptn(ji,jj,ik) * tmask(ji,jj,1)  ! masked T and S at ocean bottom
761            zsnb(ji,jj) = psn(ji,jj,ik) * tmask(ji,jj,1)
762            zdep(ji,jj) = fsdept(ji,jj,ik)                ! depth of the ocean bottom T-level
763#  if ! defined key_vectopt_loop
764         END DO
765#  endif
766      END DO
767
768      !!----------------------------------------------------------------------
769      ! 1. Criteria of additional bottom diffusivity: grad(rho).grad(h)<0
770      ! --------------------------------------------
771      ! Sign of the local density gradient along the i- and j-slopes
772      ! multiplied by the slope of the ocean bottom
773
774      SELECT CASE ( neos )
775
776      CASE ( 0 )                 ! Jackett and McDougall (1994) formulation
777
778#  if defined key_vectopt_loop
779      jj = 1
780      DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
781#  else
782      DO jj = 1, jpjm1
783         DO ji = 1, jpim1
784#  endif
785            ! temperature, salinity anomalie and depth
786            ztx = 0.5 * ( ztnb(ji,jj) + ztnb(ji+1,jj) )
787            zsx = 0.5 * ( zsnb(ji,jj) + zsnb(ji+1,jj) ) - 35.0
788            zhx = 0.5 * ( zdep(ji,jj) + zdep(ji+1,jj) )
789            !
790            zty = 0.5 * ( ztnb(ji,jj+1) + ztnb(ji,jj) )
791            zsy = 0.5 * ( zsnb(ji,jj+1) + zsnb(ji,jj) ) - 35.0
792            zhy = 0.5 * ( zdep(ji,jj+1) + zdep(ji,jj) )
793            ! masked ratio alpha/beta
794            zalbetx = fsalbt( ztx, zsx, zhx ) * umask(ji,jj,1)
795            zalbety = fsalbt( zty, zsy, zhy ) * vmask(ji,jj,1)
796            ! local density gradient along i-bathymetric slope
797            zgdrhox = zalbetx * ( ztnb(ji+1,jj) - ztnb(ji,jj) )   &
798                   -            ( zsnb(ji+1,jj) - zsnb(ji,jj) )
799            ! local density gradient along j-bathymetric slope
800            zgdrhoy = zalbety * ( ztnb(ji,jj+1) - ztnb(ji,jj) )   &
801                   -            ( zsnb(ji,jj+1) - zsnb(ji,jj) )
802            ! sign of local i-gradient of density multiplied by the i-slope
803            pbblx(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhox * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
804            ! sign of local j-gradient of density multiplied by the j-slope
805            pbbly(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhoy * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
806#  if ! defined key_vectopt_loop
807         END DO
808#  endif
809      END DO
810
811      CASE ( 1 )               ! Linear formulation function of temperature only
812                               !
813#  if defined key_vectopt_loop
814      jj = 1
815      DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
816#  else
817      DO jj = 1, jpjm1
818         DO ji = 1, jpim1
819#  endif
820            ! local 'density/temperature' gradient along i-bathymetric slope
821            zgdrhox =  ztnb(ji+1,jj) - ztnb(ji,jj)
822            ! local density gradient along j-bathymetric slope
823            zgdrhoy =  ztnb(ji,jj+1) - ztnb(ji,jj)
824            ! sign of local i-gradient of density multiplied by the i-slope
825            pbblx(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhox * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
826            ! sign of local j-gradient of density multiplied by the j-slope
827            pbbly(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhoy * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
828#  if ! defined key_vectopt_loop
829         END DO
830#  endif
831      END DO
832
833      CASE ( 2 )               ! Linear formulation function of temperature and salinity
834
835#  if defined key_vectopt_loop
836      jj = 1
837      DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
838#  else
839      DO jj = 1, jpjm1
840         DO ji = 1, jpim1
841#  endif     
842            ! local density gradient along i-bathymetric slope
843            zgdrhox = - ( rbeta*( zsnb(ji+1,jj) - zsnb(ji,jj) )   &
844                      -  ralpha*( ztnb(ji+1,jj) - ztnb(ji,jj) ) )
845            ! local density gradient along j-bathymetric slope
846            zgdrhoy = - ( rbeta*( zsnb(ji,jj+1) - zsnb(ji,jj) )   &
847                      -  ralpha*( ztnb(ji,jj+1) - ztnb(ji,jj) ) )
848            ! sign of local i-gradient of density multiplied by the i-slope
849            pbblx(ji,jj) = 0.5 - SIGN( 0.5, - zgdrhox * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
850            ! sign of local j-gradient of density multiplied by the j-slope
851            pbbly(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhoy * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
852#  if ! defined key_vectopt_loop
853         END DO
854#  endif
855      END DO
856
857      CASE DEFAULT
858
859         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
860         CALL ctl_stop(ctmp1)
861
862      END SELECT
863   
864      ! Lateral boundary conditions
865      CALL lbc_lnk( pbblx, 'U', 1. )
866      CALL lbc_lnk( pbbly, 'V', 1. )
867
868   END SUBROUTINE bbl_sign
869
870#endif
871
872   SUBROUTINE swap_dyn_data
873      !!----------------------------------------------------------------------
874      !!                    ***  ROUTINE swap_dyn_data  ***
875      !!
876      !! ** Purpose :   swap array data
877      !!
878      !! History :
879      !!   9.0  !  07-09  (C. Ethe)
880      !!----------------------------------------------------------------------
881
882
883      ! swap from record 2 to 1
884      tdta   (:,:,:,1) = tdta   (:,:,:,2)
885      sdta   (:,:,:,1) = sdta   (:,:,:,2)
886      avtdta (:,:,:,1) = avtdta (:,:,:,2)
887      udta   (:,:,:,1) = udta   (:,:,:,2)
888      vdta   (:,:,:,1) = vdta   (:,:,:,2)
889      wdta   (:,:,:,1) = wdta   (:,:,:,2)
890#if defined key_trc_diatrd
891      hdivdta(:,:,:,1) = hdivdta(:,:,:,2)
892#endif
893
894#if defined key_ldfslp
895      uslpdta (:,:,:,1) = uslpdta (:,:,:,2)
896      vslpdta (:,:,:,1) = vslpdta (:,:,:,2)
897      wslpidta(:,:,:,1) = wslpidta(:,:,:,2)
898      wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2)
899#endif
900      hmlddta(:,:,1) = hmlddta(:,:,2) 
901      wspddta(:,:,1) = wspddta(:,:,2) 
902      frlddta(:,:,1) = frlddta(:,:,2) 
903      empdta (:,:,1) = empdta (:,:,2) 
904      qsrdta (:,:,1) = qsrdta (:,:,2) 
905
906#if ! defined key_off_degrad && defined key_traldf_c2d
907      ahtwdta(:,:,1) = ahtwdta(:,:,2)
908#  if defined key_trcldf_eiv
909      aeiwdta(:,:,1) = aeiwdta(:,:,2)
910#  endif
911#endif
912
913#if defined key_off_degrad
914      ahtudta(:,:,:,1) = ahtudta(:,:,:,2)
915      ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2)
916      ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2)
917#  if defined key_trcldf_eiv
918      aeiudta(:,:,:,1) = aeiudta(:,:,:,2)
919      aeivdta(:,:,:,1) = aeivdta(:,:,:,2)
920      aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2)
921#  endif
922#endif
923
924#if defined key_trcbbl_dif || defined key_trcbbl_adv
925      bblxdta(:,:,1) = bblxdta(:,:,2)
926      bblydta(:,:,1) = bblydta(:,:,2)
927#endif
928
929   END SUBROUTINE swap_dyn_data
930
931   SUBROUTINE assign_dyn_data
932      !!----------------------------------------------------------------------
933      !!                    ***  ROUTINE assign_dyn_data  ***
934      !!
935      !! ** Purpose :   Assign dynamical data to the data that have been read
936      !!                without time interpolation
937      !!
938      !!----------------------------------------------------------------------
939     
940      tn (:,:,:) = tdta  (:,:,:,2)
941      sn (:,:,:) = sdta  (:,:,:,2)
942      avt(:,:,:) = avtdta(:,:,:,2)
943     
944      un (:,:,:) = udta  (:,:,:,2) 
945      vn (:,:,:) = vdta  (:,:,:,2)
946      wn (:,:,:) = wdta  (:,:,:,2)
947
948#if defined key_trc_diatrd
949      hdivn(:,:,:) = hdivdta(:,:,:,2)
950#endif
951
952#if defined key_zdfddm
953      avs(:,:,:)   = avtdta (:,:,:,2)
954#endif
955
956     
957#if defined key_ldfslp
958      uslp (:,:,:) = uslpdta (:,:,:,2) 
959      vslp (:,:,:) = vslpdta (:,:,:,2) 
960      wslpi(:,:,:) = wslpidta(:,:,:,2) 
961      wslpj(:,:,:) = wslpjdta(:,:,:,2) 
962#endif
963
964      hmld(:,:) = hmlddta(:,:,2) 
965      wndm(:,:) = wspddta(:,:,2) 
966      fr_i(:,:) = frlddta(:,:,2) 
967      emp (:,:) = empdta (:,:,2) 
968      emps(:,:) = emp(:,:) 
969      qsr (:,:) = qsrdta (:,:,2) 
970
971#if ! defined key_off_degrad && defined key_traldf_c2d   
972      ahtw(:,:) = ahtwdta(:,:,2)
973#  if defined key_trcldf_eiv
974      aeiw(:,:) = aeiwdta(:,:,2)
975#  endif
976#endif
977     
978#if defined key_off_degrad
979      ahtu(:,:,:) = ahtudta(:,:,:,2)
980      ahtv(:,:,:) = ahtvdta(:,:,:,2)
981      ahtw(:,:,:) = ahtwdta(:,:,:,2)
982#  if defined key_trcldf_eiv
983      aeiu(:,:,:) = aeiudta(:,:,:,2)
984      aeiv(:,:,:) = aeivdta(:,:,:,2)
985      aeiw(:,:,:) = aeiwdta(:,:,:,2)
986#  endif
987     
988#endif
989     
990#if defined key_trcbbl_dif ||  defined key_trcbbl_adv
991      bblx(:,:) = bblxdta(:,:,2)
992      bbly(:,:) = bblydta(:,:,2)
993#endif
994
995   END SUBROUTINE assign_dyn_data
996
997   SUBROUTINE linear_interp_dyn_data( pweigh )
998      !!----------------------------------------------------------------------
999      !!                    ***  ROUTINE linear_interp_dyn_data  ***
1000      !!
1001      !! ** Purpose :   linear interpolation of data
1002      !!
1003      !!----------------------------------------------------------------------
1004      !! * Argument
1005      REAL(wp), INTENT( in ) ::   pweigh       ! weigh
1006
1007      !! * Local declarations
1008      REAL(wp) :: zweighm1
1009      !!----------------------------------------------------------------------
1010
1011      zweighm1 = 1. - pweigh
1012     
1013      tn (:,:,:) = zweighm1 * tdta  (:,:,:,1) + pweigh * tdta  (:,:,:,2)
1014      sn (:,:,:) = zweighm1 * sdta  (:,:,:,1) + pweigh * sdta  (:,:,:,2)
1015      avt(:,:,:) = zweighm1 * avtdta(:,:,:,1) + pweigh * avtdta(:,:,:,2)
1016     
1017      un (:,:,:) = zweighm1 * udta  (:,:,:,1) + pweigh * udta  (:,:,:,2) 
1018      vn (:,:,:) = zweighm1 * vdta  (:,:,:,1) + pweigh * vdta  (:,:,:,2)
1019      wn (:,:,:) = zweighm1 * wdta  (:,:,:,1) + pweigh * wdta  (:,:,:,2)
1020
1021#if defined key_trc_diatrd
1022      hdivn(:,:,:) = zweighm1 * hdivdta(:,:,:,1) + pweigh * hdivdta(:,:,:,2)
1023#endif
1024
1025#if defined key_zdfddm
1026      avs(:,:,:)   = zweighm1 * avtdta (:,:,:,1) + pweigh * avtdta (:,:,:,2)
1027#endif
1028
1029     
1030#if defined key_ldfslp
1031      uslp (:,:,:) = zweighm1 * uslpdta (:,:,:,1) + pweigh * uslpdta (:,:,:,2) 
1032      vslp (:,:,:) = zweighm1 * vslpdta (:,:,:,1) + pweigh * vslpdta (:,:,:,2) 
1033      wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + pweigh * wslpidta(:,:,:,2) 
1034      wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + pweigh * wslpjdta(:,:,:,2) 
1035#endif
1036
1037      hmld(:,:) = zweighm1 * hmlddta(:,:,1) + pweigh  * hmlddta(:,:,2) 
1038      wndm(:,:) = zweighm1 * wspddta(:,:,1) + pweigh  * wspddta(:,:,2) 
1039      fr_i(:,:) = zweighm1 * frlddta(:,:,1) + pweigh  * frlddta(:,:,2) 
1040      emp (:,:) = zweighm1 * empdta (:,:,1) + pweigh  * empdta (:,:,2) 
1041      emps(:,:) = emp(:,:) 
1042      qsr (:,:) = zweighm1 * qsrdta (:,:,1) + pweigh  * qsrdta (:,:,2) 
1043
1044#if ! defined key_off_degrad && defined key_traldf_c2d   
1045      ahtw(:,:) = zweighm1 * ahtwdta(:,:,1) + pweigh * ahtwdta(:,:,2)
1046#  if defined key_trcldf_eiv
1047      aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + pweigh * aeiwdta(:,:,2)
1048#  endif
1049#endif
1050     
1051#if defined key_off_degrad
1052      ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + pweigh * ahtudta(:,:,:,2)
1053      ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + pweigh * ahtvdta(:,:,:,2)
1054      ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + pweigh * ahtwdta(:,:,:,2)
1055#  if defined key_trcldf_eiv
1056      aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + pweigh * aeiudta(:,:,:,2)
1057      aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + pweigh * aeivdta(:,:,:,2)
1058      aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + pweigh * aeiwdta(:,:,:,2)
1059#  endif
1060#endif
1061     
1062#if defined key_trcbbl_dif   ||   defined key_trcbbl_adv
1063      bblx(:,:) = zweighm1 * bblxdta(:,:,1) + pweigh * bblxdta(:,:,2)
1064      bbly(:,:) = zweighm1 * bblydta(:,:,1) + pweigh * bblydta(:,:,2)
1065#endif
1066
1067   END SUBROUTINE linear_interp_dyn_data
1068
1069END MODULE dtadyn
Note: See TracBrowser for help on using the repository browser.