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.
ldftra.F90 in branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/LDF – NEMO

source: branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/LDF/ldftra.F90 @ 6748

Last change on this file since 6748 was 6748, checked in by mocavero, 8 years ago

GYRE hybrid parallelization

  • Property svn:keywords set to Id
File size: 41.8 KB
Line 
1MODULE ldftra
2   !!======================================================================
3   !!                       ***  MODULE  ldftra  ***
4   !! Ocean physics:  lateral diffusivity coefficients
5   !!=====================================================================
6   !! History :       ! 1997-07  (G. Madec)  from inimix.F split in 2 routines
7   !!   NEMO     1.0  ! 2002-09  (G. Madec)  F90: Free form and module
8   !!            2.0  ! 2005-11  (G. Madec) 
9   !!            3.7  ! 2013-12  (F. Lemarie, G. Madec)  restructuration/simplification of aht/aeiv specification,
10   !!                 !                                  add velocity dependent coefficient and optional read in file
11   !!----------------------------------------------------------------------
12
13   !!----------------------------------------------------------------------
14   !!   ldf_tra_init : initialization, namelist read, and parameters control
15   !!   ldf_tra      : update lateral eddy diffusivity coefficients at each time step
16   !!   ldf_eiv_init : initialization of the eiv coeff. from namelist choices
17   !!   ldf_eiv      : time evolution of the eiv coefficients (function of the growth rate of baroclinic instability)
18   !!   ldf_eiv_trp  : add to the input ocean transport the contribution of the EIV parametrization
19   !!   ldf_eiv_dia  : diagnose the eddy induced velocity from the eiv streamfunction
20   !!----------------------------------------------------------------------
21   USE oce             ! ocean dynamics and tracers
22   USE dom_oce         ! ocean space and time domain
23   USE phycst          ! physical constants
24   USE ldfslp          ! lateral diffusion: slope of iso-neutral surfaces
25   USE ldfc1d_c2d      ! lateral diffusion: 1D & 2D cases
26   USE diaar5, ONLY:   lk_diaar5
27   !
28   USE trc_oce, ONLY: lk_offline ! offline flag
29   USE in_out_manager  ! I/O manager
30   USE iom             ! I/O module for ehanced bottom friction file
31   USE lib_mpp         ! distribued memory computing library
32   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
33   USE wrk_nemo        ! work arrays
34   USE timing          ! timing
35
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   ldf_tra_init   ! called by nemogcm.F90
40   PUBLIC   ldf_tra        ! called by step.F90
41   PUBLIC   ldf_eiv_init   ! called by nemogcm.F90
42   PUBLIC   ldf_eiv        ! called by step.F90
43   PUBLIC   ldf_eiv_trp    ! called by traadv.F90
44   PUBLIC   ldf_eiv_dia    ! called by traldf_iso and traldf_iso_triad.F90
45   
46   !                                   !!* Namelist namtra_ldf : lateral mixing on tracers *
47   !                                    != Operator type =!
48   LOGICAL , PUBLIC ::   ln_traldf_lap       !: laplacian operator
49   LOGICAL , PUBLIC ::   ln_traldf_blp       !: bilaplacian operator
50   !                                    != Direction of action =!
51   LOGICAL , PUBLIC ::   ln_traldf_lev       !: iso-level direction
52   LOGICAL , PUBLIC ::   ln_traldf_hor       !: horizontal (geopotential) direction
53!  LOGICAL , PUBLIC ::   ln_traldf_iso       !: iso-neutral direction                    (see ldfslp)
54!  LOGICAL , PUBLIC ::   ln_traldf_triad     !: griffies triad scheme                    (see ldfslp)
55   LOGICAL , PUBLIC ::   ln_traldf_msc       !: Method of Stabilizing Correction
56!  LOGICAL , PUBLIC ::   ln_triad_iso        !: pure horizontal mixing in ML             (see ldfslp)
57!  LOGICAL , PUBLIC ::   ln_botmix_triad     !: mixing on bottom                         (see ldfslp)
58!  REAL(wp), PUBLIC ::   rn_sw_triad         !: =1/0 switching triad / all 4 triads used (see ldfslp)
59!  REAL(wp), PUBLIC ::   rn_slpmax           !: slope limit                              (see ldfslp)
60   !                                    !=  Coefficients =!
61   INTEGER , PUBLIC ::   nn_aht_ijk_t        !: choice of time & space variations of the lateral eddy diffusivity coef.
62   REAL(wp), PUBLIC ::   rn_aht_0            !:   laplacian lateral eddy diffusivity [m2/s]
63   REAL(wp), PUBLIC ::   rn_bht_0            !: bilaplacian lateral eddy diffusivity [m4/s]
64
65   !                                   !!* Namelist namtra_ldfeiv : eddy induced velocity param. *
66   !                                    != Use/diagnose eiv =!
67   LOGICAL , PUBLIC ::   ln_ldfeiv           !: eddy induced velocity flag
68   LOGICAL , PUBLIC ::   ln_ldfeiv_dia       !: diagnose & output eiv streamfunction and velocity (IOM)
69   !                                    != Coefficients =!
70   INTEGER , PUBLIC ::   nn_aei_ijk_t        !: choice of time/space variation of the eiv coeff.
71   REAL(wp), PUBLIC ::   rn_aeiv_0           !: eddy induced velocity coefficient [m2/s]
72   
73   LOGICAL , PUBLIC ::   l_ldftra_time = .FALSE.   !: flag for time variation of the lateral eddy diffusivity coef.
74   LOGICAL , PUBLIC ::   l_ldfeiv_time = .FALSE.   ! flag for time variation of the eiv coef.
75
76   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ahtu, ahtv   !: eddy diffusivity coef. at U- and V-points   [m2/s]
77   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   aeiu, aeiv   !: eddy induced velocity coeff.                [m2/s]
78
79   REAL(wp) ::   r1_4  = 0.25_wp          ! =1/4
80   REAL(wp) ::   r1_12 = 1._wp / 12._wp   ! =1/12
81
82   !! * Substitutions
83#  include "vectopt_loop_substitute.h90"
84   !!----------------------------------------------------------------------
85   !! NEMO/OPA 3.7 , NEMO Consortium (2015)
86   !! $Id$
87   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
88   !!----------------------------------------------------------------------
89CONTAINS
90
91   SUBROUTINE ldf_tra_init
92      !!----------------------------------------------------------------------
93      !!                  ***  ROUTINE ldf_tra_init  ***
94      !!
95      !! ** Purpose :   initializations of the tracer lateral mixing coeff.
96      !!
97      !! ** Method  : * the eddy diffusivity coef. specification depends on:
98      !!
99      !!    ln_traldf_lap = T     laplacian operator
100      !!    ln_traldf_blp = T   bilaplacian operator
101      !!
102      !!    nn_aht_ijk_t  =  0 => = constant
103      !!                  !
104      !!                  = 10 => = F(z) : constant with a reduction of 1/4 with depth
105      !!                  !
106      !!                  =-20 => = F(i,j)   = shape read in 'eddy_diffusivity.nc' file
107      !!                  = 20    = F(i,j)   = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
108      !!                  = 21    = F(i,j,t) = F(growth rate of baroclinic instability)
109      !!                  !
110      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_diffusivity.nc' file
111      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
112      !!                  = 31    = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
113      !!                                                          or |u|e^3/12 bilaplacian operator )
114      !!              * initialisation of the eddy induced velocity coefficient by a call to ldf_eiv_init
115      !!           
116      !! ** action  : ahtu, ahtv initialized once for all or l_ldftra_time set to true
117      !!              aeiu, aeiv initialized once for all or l_ldfeiv_time set to true
118      !!----------------------------------------------------------------------
119      INTEGER  ::   jk                ! dummy loop indices
120      INTEGER  ::   ierr, inum, ios   ! local integer
121      REAL(wp) ::   zah0              ! local scalar
122      !
123      NAMELIST/namtra_ldf/ ln_traldf_lap, ln_traldf_blp  ,                   &   ! type of operator
124         &                 ln_traldf_lev, ln_traldf_hor  , ln_traldf_triad,  &   ! acting direction of the operator
125         &                 ln_traldf_iso, ln_traldf_msc  ,  rn_slpmax     ,  &   ! option for iso-neutral operator
126         &                 ln_triad_iso , ln_botmix_triad, rn_sw_triad    ,  &   ! option for triad operator
127         &                 rn_aht_0     , rn_bht_0       , nn_aht_ijk_t          ! lateral eddy coefficient
128      !!----------------------------------------------------------------------
129      !
130      !  Choice of lateral tracer physics
131      ! =================================
132      !
133      REWIND( numnam_ref )              ! Namelist namtra_ldf in reference namelist : Lateral physics on tracers
134      READ  ( numnam_ref, namtra_ldf, IOSTAT = ios, ERR = 901)
135901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_ldf in reference namelist', lwp )
136      !
137      REWIND( numnam_cfg )              ! Namelist namtra_ldf in configuration namelist : Lateral physics on tracers
138      READ  ( numnam_cfg, namtra_ldf, IOSTAT = ios, ERR = 902 )
139902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_ldf in configuration namelist', lwp )
140      IF(lwm) WRITE ( numond, namtra_ldf )
141      !
142      IF(lwp) THEN                      ! control print
143         WRITE(numout,*)
144         WRITE(numout,*) 'ldf_tra_init : lateral tracer physics'
145         WRITE(numout,*) '~~~~~~~~~~~~ '
146         WRITE(numout,*) '   Namelist namtra_ldf : lateral mixing parameters (type, direction, coefficients)'
147         !
148         WRITE(numout,*) '      type :'
149         WRITE(numout,*) '         laplacian operator                      ln_traldf_lap   = ', ln_traldf_lap
150         WRITE(numout,*) '         bilaplacian operator                    ln_traldf_blp   = ', ln_traldf_blp
151         !
152         WRITE(numout,*) '      direction of action :'
153         WRITE(numout,*) '         iso-level                               ln_traldf_lev   = ', ln_traldf_lev
154         WRITE(numout,*) '         horizontal (geopotential)               ln_traldf_hor   = ', ln_traldf_hor
155         WRITE(numout,*) '         iso-neutral Madec operator              ln_traldf_iso   = ', ln_traldf_iso
156         WRITE(numout,*) '         iso-neutral triad operator              ln_traldf_triad = ', ln_traldf_triad
157         WRITE(numout,*) '            iso-neutral (Method of Stab. Corr.)  ln_traldf_msc   = ', ln_traldf_msc
158         WRITE(numout,*) '            maximum isoppycnal slope             rn_slpmax       = ', rn_slpmax
159         WRITE(numout,*) '            pure lateral mixing in ML            ln_triad_iso    = ', ln_triad_iso
160         WRITE(numout,*) '            switching triad or not               rn_sw_triad     = ', rn_sw_triad
161         WRITE(numout,*) '            lateral mixing on bottom             ln_botmix_triad = ', ln_botmix_triad
162         !
163         WRITE(numout,*) '      coefficients :'
164         WRITE(numout,*) '         lateral eddy diffusivity   (lap case)   rn_aht_0        = ', rn_aht_0
165         WRITE(numout,*) '         lateral eddy diffusivity (bilap case)   rn_bht_0        = ', rn_bht_0
166         WRITE(numout,*) '         type of time-space variation            nn_aht_ijk_t    = ', nn_aht_ijk_t
167      ENDIF
168      !
169      !                                ! Parameter control
170      !
171      IF( .NOT.ln_traldf_lap .AND. .NOT.ln_traldf_blp ) THEN
172         IF(lwp) WRITE(numout,*) '   No diffusive operator selected. ahtu and ahtv are not allocated'
173         l_ldftra_time = .FALSE.
174         RETURN
175      ENDIF
176      !
177      IF( ln_traldf_blp .AND. ( ln_traldf_iso .OR. ln_traldf_triad) ) THEN     ! iso-neutral bilaplacian need MSC
178         IF( .NOT.ln_traldf_msc )   CALL ctl_stop( 'tra_ldf_init: iso-neutral bilaplacian requires ln_traldf_msc=.true.' )
179      ENDIF
180      !
181      !  Space/time variation of eddy coefficients
182      ! ===========================================
183      !                                               ! allocate the aht arrays
184      ALLOCATE( ahtu(jpi,jpj,jpk) , ahtv(jpi,jpj,jpk) , STAT=ierr )
185      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_tra_init: failed to allocate arrays')
186      !
187!$OMP PARALLEL WORKSHARE
188      ahtu(:,:,jpk) = 0._wp                           ! last level always 0 
189      ahtv(:,:,jpk) = 0._wp
190!$OMP END PARALLEL WORKSHARE
191      !
192      !                                               ! value of eddy mixing coef.
193      IF    ( ln_traldf_lap ) THEN   ;   zah0 =      rn_aht_0        !   laplacian operator
194      ELSEIF( ln_traldf_blp ) THEN   ;   zah0 = ABS( rn_bht_0 )      ! bilaplacian operator
195      ENDIF
196      !
197      l_ldftra_time = .FALSE.                         ! no time variation except in case defined below
198      !
199      IF( ln_traldf_lap .OR. ln_traldf_blp ) THEN     ! only if a lateral diffusion operator is used
200         !
201         SELECT CASE(  nn_aht_ijk_t  )                   ! Specification of space time variations of ehtu, ahtv
202         !
203         CASE(   0  )      !==  constant  ==!
204            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = constant = ', rn_aht_0
205!$OMP PARALLEL WORKSHARE
206            ahtu(:,:,:) = zah0 * umask(:,:,:)
207            ahtv(:,:,:) = zah0 * vmask(:,:,:)
208!$OMP END PARALLEL WORKSHARE
209            !
210         CASE(  10  )      !==  fixed profile  ==!
211            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( depth )'
212!$OMP PARALLEL WORKSHARE
213            ahtu(:,:,1) = zah0 * umask(:,:,1)                      ! constant surface value
214            ahtv(:,:,1) = zah0 * vmask(:,:,1)
215!$OMP END PARALLEL WORKSHARE
216            CALL ldf_c1d( 'TRA', r1_4, ahtu(:,:,1), ahtv(:,:,1), ahtu, ahtv )
217            !
218         CASE ( -20 )      !== fixed horizontal shape read in file  ==!
219            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F(i,j) read in eddy_diffusivity.nc file'
220            CALL iom_open( 'eddy_diffusivity_2D.nc', inum )
221            CALL iom_get ( inum, jpdom_data, 'ahtu_2D', ahtu(:,:,1) )
222            CALL iom_get ( inum, jpdom_data, 'ahtv_2D', ahtv(:,:,1) )
223            CALL iom_close( inum )
224!$OMP PARALLEL DO schedule(static) private(jk)
225            DO jk = 2, jpkm1
226               ahtu(:,:,jk) = ahtu(:,:,1) * umask(:,:,jk)
227               ahtv(:,:,jk) = ahtv(:,:,1) * vmask(:,:,jk)
228            END DO
229            !
230         CASE(  20  )      !== fixed horizontal shape  ==!
231            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( e1, e2 ) or F( e1^3, e2^3 ) (lap or blp case)'
232            IF( ln_traldf_lap )   CALL ldf_c2d( 'TRA', 'LAP', zah0, ahtu, ahtv )    ! surface value proportional to scale factor
233            IF( ln_traldf_blp )   CALL ldf_c2d( 'TRA', 'BLP', zah0, ahtu, ahtv )    ! surface value proportional to scale factor
234            !
235         CASE(  21  )      !==  time varying 2D field  ==!
236            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( latitude, longitude, time )'
237            IF(lwp) WRITE(numout,*) '                              = F( growth rate of baroclinic instability )'
238            IF(lwp) WRITE(numout,*) '                              min value = 0.1 * rn_aht_0'
239            IF(lwp) WRITE(numout,*) '                              max value = rn_aht_0 (rn_aeiv_0 if nn_aei_ijk_t=21)'
240            IF(lwp) WRITE(numout,*) '                              increased to rn_aht_0 within 20N-20S'
241            !
242            l_ldftra_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
243            !
244            IF( ln_traldf_blp ) THEN
245               CALL ctl_stop( 'ldf_tra_init: aht=F(growth rate of baroc. insta.) incompatible with bilaplacian operator' )
246            ENDIF
247            !
248         CASE( -30  )      !== fixed 3D shape read in file  ==!
249            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F(i,j,k) read in eddy_diffusivity.nc file'
250            CALL iom_open( 'eddy_diffusivity_3D.nc', inum )
251            CALL iom_get ( inum, jpdom_data, 'ahtu_3D', ahtu )
252            CALL iom_get ( inum, jpdom_data, 'ahtv_3D', ahtv )
253            CALL iom_close( inum )
254!$OMP PARALLEL DO schedule(static) private(jk)
255            DO jk = 1, jpkm1
256               ahtu(:,:,jk) = ahtu(:,:,jk) * umask(:,:,jk)
257               ahtv(:,:,jk) = ahtv(:,:,jk) * vmask(:,:,jk)
258            END DO
259            !
260         CASE(  30  )      !==  fixed 3D shape  ==!
261            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( latitude, longitude, depth )'
262            IF( ln_traldf_lap )   CALL ldf_c2d( 'TRA', 'LAP', zah0, ahtu, ahtv )    ! surface value proportional to scale factor
263            IF( ln_traldf_blp )   CALL ldf_c2d( 'TRA', 'BLP', zah0, ahtu, ahtv )    ! surface value proportional to scale factor
264            !                                                    ! reduction with depth
265            CALL ldf_c1d( 'TRA', r1_4, ahtu(:,:,1), ahtv(:,:,1), ahtu, ahtv )
266            !
267         CASE(  31  )      !==  time varying 3D field  ==!
268            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( latitude, longitude, depth , time )'
269            IF(lwp) WRITE(numout,*) '                                proportional to the velocity : |u|e/12 or |u|e^3/12'
270            !
271            l_ldftra_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
272            !
273         CASE DEFAULT
274            CALL ctl_stop('ldf_tra_init: wrong choice for nn_aht_ijk_t, the type of space-time variation of aht')
275         END SELECT
276         !
277         IF( ln_traldf_blp .AND. .NOT. l_ldftra_time ) THEN
278!$OMP PARALLEL WORKSHARE
279            ahtu(:,:,:) = SQRT( ahtu(:,:,:) )
280            ahtv(:,:,:) = SQRT( ahtv(:,:,:) )
281!$OMP END PARALLEL WORKSHARE
282         ENDIF
283         !
284      ENDIF
285      !
286   END SUBROUTINE ldf_tra_init
287
288
289   SUBROUTINE ldf_tra( kt )
290      !!----------------------------------------------------------------------
291      !!                  ***  ROUTINE ldf_tra  ***
292      !!
293      !! ** Purpose :   update at kt the tracer lateral mixing coeff. (aht and aeiv)
294      !!
295      !! ** Method  :   time varying eddy diffusivity coefficients:
296      !!
297      !!    nn_aei_ijk_t = 21    aeiu, aeiv = F(i,j,  t) = F(growth rate of baroclinic instability)
298      !!                                                   with a reduction to 0 in vicinity of the Equator
299      !!    nn_aht_ijk_t = 21    ahtu, ahtv = F(i,j,  t) = F(growth rate of baroclinic instability)
300      !!
301      !!                 = 31    ahtu, ahtv = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
302      !!                                                                     or |u|e^3/12 bilaplacian operator )
303      !!
304      !! ** action  :   ahtu, ahtv   update at each time step   
305      !!                aeiu, aeiv      -       -     -    -   (if ln_ldfeiv=T)
306      !!----------------------------------------------------------------------
307      INTEGER, INTENT(in) ::   kt   ! time step
308      !
309      INTEGER  ::   ji, jj, jk   ! dummy loop indices
310      REAL(wp) ::   zaht, zaht_min, z1_f20       ! local scalar
311      !!----------------------------------------------------------------------
312      !
313      IF( nn_aei_ijk_t == 21 ) THEN       ! eddy induced velocity coefficients
314         !                                ! =F(growth rate of baroclinic instability)
315         !                                ! max value rn_aeiv_0 ; decreased to 0 within 20N-20S
316         CALL ldf_eiv( kt, rn_aeiv_0, aeiu, aeiv )
317         IF(lwp .AND. kt<=nit000+20 )   WRITE(numout,*) ' kt , ldf_eiv appel', kt
318      ENDIF
319      !
320      SELECT CASE(  nn_aht_ijk_t  )       ! Eddy diffusivity coefficients
321      !
322      CASE(  21  )       !==  time varying 2D field  ==!   = F( growth rate of baroclinic instability )
323         !                                             !   min value rn_aht_0 / 10
324         !                                             !   max value rn_aht_0 (rn_aeiv_0 if nn_aei_ijk_t=21)
325         !                                             !   increase to rn_aht_0 within 20N-20S
326         IF( nn_aei_ijk_t /= 21 ) THEN
327            CALL ldf_eiv( kt, rn_aht_0, ahtu, ahtv )
328            IF(lwp .AND. kt<=nit000+20 )   WRITE(numout,*) ' kt , ldf_eiv appel  2', kt
329         ELSE
330            ahtu(:,:,1) = aeiu(:,:,1)
331            ahtv(:,:,1) = aeiv(:,:,1)
332            IF(lwp .AND. kt<=nit000+20 )   WRITE(numout,*) ' kt , ahtu=aeiu', kt
333         ENDIF
334         !
335         z1_f20   = 1._wp / (  2._wp * omega * SIN( rad * 20._wp )  )      ! 1 / ff(20 degrees)   
336         zaht_min = 0.2_wp * rn_aht_0                                      ! minimum value for aht
337         DO jj = 1, jpj
338            DO ji = 1, jpi
339               zaht = ( 1._wp -  MIN( 1._wp , ABS( ff(ji,jj) * z1_f20 ) ) ) * ( rn_aht_0 - zaht_min )
340               ahtu(ji,jj,1) = (  MAX( zaht_min, ahtu(ji,jj,1) ) + zaht  ) * umask(ji,jj,1)     ! min value zaht_min
341               ahtv(ji,jj,1) = (  MAX( zaht_min, ahtv(ji,jj,1) ) + zaht  ) * vmask(ji,jj,1)     ! increase within 20S-20N
342            END DO
343         END DO
344         DO jk = 2, jpkm1                             ! deeper value = surface value
345            ahtu(:,:,jk) = ahtu(:,:,1) * umask(:,:,jk)
346            ahtv(:,:,jk) = ahtv(:,:,1) * vmask(:,:,jk)
347         END DO
348         !
349      CASE(  31  )       !==  time varying 3D field  ==!   = F( local velocity )
350         IF( ln_traldf_lap     ) THEN          !   laplacian operator |u| e /12
351            DO jk = 1, jpkm1
352               ahtu(:,:,jk) = ABS( ub(:,:,jk) ) * e1u(:,:) * r1_12
353               ahtv(:,:,jk) = ABS( vb(:,:,jk) ) * e2v(:,:) * r1_12
354            END DO
355         ELSEIF( ln_traldf_blp ) THEN      ! bilaplacian operator      sqrt( |u| e^3 /12 ) = sqrt( |u| e /12 ) * e
356            DO jk = 1, jpkm1
357               ahtu(:,:,jk) = SQRT(  ABS( ub(:,:,jk) ) * e1u(:,:) * r1_12  ) * e1u(:,:)
358               ahtv(:,:,jk) = SQRT(  ABS( vb(:,:,jk) ) * e2v(:,:) * r1_12  ) * e2v(:,:)
359            END DO
360         ENDIF
361         !
362      END SELECT
363      !
364      IF( .NOT.lk_offline ) THEN
365         CALL iom_put( "ahtu_2d", ahtu(:,:,1) )   ! surface u-eddy diffusivity coeff.
366         CALL iom_put( "ahtv_2d", ahtv(:,:,1) )   ! surface v-eddy diffusivity coeff.
367         CALL iom_put( "ahtu_3d", ahtu(:,:,:) )   ! 3D      u-eddy diffusivity coeff.
368         CALL iom_put( "ahtv_3d", ahtv(:,:,:) )   ! 3D      v-eddy diffusivity coeff.
369         !
370!!gm  : THE IF below is to be checked (comes from Seb)
371         IF( ln_ldfeiv ) THEN
372           CALL iom_put( "aeiu_2d", aeiu(:,:,1) )   ! surface u-EIV coeff.
373           CALL iom_put( "aeiv_2d", aeiv(:,:,1) )   ! surface v-EIV coeff.
374           CALL iom_put( "aeiu_3d", aeiu(:,:,:) )   ! 3D      u-EIV coeff.
375           CALL iom_put( "aeiv_3d", aeiv(:,:,:) )   ! 3D      v-EIV coeff.
376         ENDIF     
377      ENDIF
378      !
379   END SUBROUTINE ldf_tra
380
381
382   SUBROUTINE ldf_eiv_init
383      !!----------------------------------------------------------------------
384      !!                  ***  ROUTINE ldf_eiv_init  ***
385      !!
386      !! ** Purpose :   initialization of the eiv coeff. from namelist choices.
387      !!
388      !! ** Method :
389      !!
390      !! ** Action :   aeiu , aeiv   : EIV coeff. at u- & v-points
391      !!               l_ldfeiv_time : =T if EIV coefficients vary with time
392      !!----------------------------------------------------------------------
393      INTEGER  ::   jk                ! dummy loop indices
394      INTEGER  ::   ierr, inum, ios   ! local integer
395      !
396      NAMELIST/namtra_ldfeiv/ ln_ldfeiv   , ln_ldfeiv_dia,   &    ! eddy induced velocity (eiv)
397         &                    nn_aei_ijk_t, rn_aeiv_0             ! eiv  coefficient
398      !!----------------------------------------------------------------------
399      !
400      REWIND( numnam_ref )              ! Namelist namtra_ldfeiv in reference namelist : eddy induced velocity param.
401      READ  ( numnam_ref, namtra_ldfeiv, IOSTAT = ios, ERR = 901)
402901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_ldfeiv in reference namelist', lwp )
403      !
404      REWIND( numnam_cfg )              ! Namelist namtra_ldfeiv in configuration namelist : eddy induced velocity param.
405      READ  ( numnam_cfg, namtra_ldfeiv, IOSTAT = ios, ERR = 902 )
406902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_ldfeiv in configuration namelist', lwp )
407      IF(lwm)  WRITE ( numond, namtra_ldfeiv )
408
409      IF(lwp) THEN                      ! control print
410         WRITE(numout,*)
411         WRITE(numout,*) 'ldf_eiv_init : eddy induced velocity parametrization'
412         WRITE(numout,*) '~~~~~~~~~~~~ '
413         WRITE(numout,*) '   Namelist namtra_ldfeiv : '
414         WRITE(numout,*) '      Eddy Induced Velocity (eiv) param.      ln_ldfeiv     = ', ln_ldfeiv
415         WRITE(numout,*) '      eiv streamfunction & velocity diag.     ln_ldfeiv_dia = ', ln_ldfeiv_dia
416         WRITE(numout,*) '      eddy induced velocity coef.             rn_aeiv_0     = ', rn_aeiv_0
417         WRITE(numout,*) '      type of time-space variation            nn_aei_ijk_t  = ', nn_aei_ijk_t
418         WRITE(numout,*)
419      ENDIF
420      !
421      IF( ln_ldfeiv .AND. ln_traldf_blp )   CALL ctl_stop( 'ldf_eiv_init: eddy induced velocity ONLY with laplacian diffusivity' )
422
423      !                                 ! Parameter control
424      l_ldfeiv_time = .FALSE.   
425      !
426      IF( ln_ldfeiv ) THEN                         ! allocate the aei arrays
427         ALLOCATE( aeiu(jpi,jpj,jpk), aeiv(jpi,jpj,jpk), STAT=ierr )
428         IF( ierr /= 0 )   CALL ctl_stop('STOP', 'ldf_eiv: failed to allocate arrays')
429         !
430         SELECT CASE( nn_aei_ijk_t )               ! Specification of space time variations of eaiu, aeiv
431         !
432         CASE(   0  )      !==  constant  ==!
433            IF(lwp) WRITE(numout,*) '          eddy induced velocity coef. = constant = ', rn_aeiv_0
434!$OMP PARALLEL WORKSHARE
435            aeiu(:,:,:) = rn_aeiv_0
436            aeiv(:,:,:) = rn_aeiv_0
437!$OMP END PARALLEL WORKSHARE
438            !
439         CASE(  10  )      !==  fixed profile  ==!
440            IF(lwp) WRITE(numout,*) '          eddy induced velocity coef. = F( depth )'
441!$OMP PARALLEL WORKSHARE
442            aeiu(:,:,1) = rn_aeiv_0                                ! constant surface value
443            aeiv(:,:,1) = rn_aeiv_0
444!$OMP END PARALLEL WORKSHARE
445            CALL ldf_c1d( 'TRA', r1_4, aeiu(:,:,1), aeiv(:,:,1), aeiu, aeiv )
446            !
447         CASE ( -20 )      !== fixed horizontal shape read in file  ==!
448            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F(i,j) read in eddy_diffusivity_2D.nc file'
449            CALL iom_open ( 'eddy_induced_velocity_2D.nc', inum )
450            CALL iom_get  ( inum, jpdom_data, 'aeiu', aeiu(:,:,1) )
451            CALL iom_get  ( inum, jpdom_data, 'aeiv', aeiv(:,:,1) )
452            CALL iom_close( inum )
453!$OMP PARALLEL DO schedule(static) private(jk)
454            DO jk = 2, jpk
455               aeiu(:,:,jk) = aeiu(:,:,1)
456               aeiv(:,:,jk) = aeiv(:,:,1)
457            END DO
458            !
459         CASE(  20  )      !== fixed horizontal shape  ==!
460            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( e1, e2 ) or F( e1^3, e2^3 ) (lap or bilap case)'
461            CALL ldf_c2d( 'TRA', 'LAP', rn_aeiv_0, aeiu, aeiv )    ! surface value proportional to scale factor
462            !
463         CASE(  21  )       !==  time varying 2D field  ==!
464            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( latitude, longitude, time )'
465            IF(lwp) WRITE(numout,*) '                              = F( growth rate of baroclinic instability )'
466            !
467            l_ldfeiv_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
468            !
469         CASE( -30  )      !== fixed 3D shape read in file  ==!
470            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F(i,j,k) read in eddy_diffusivity_3D.nc file'
471            CALL iom_open ( 'eddy_induced_velocity_3D.nc', inum )
472            CALL iom_get  ( inum, jpdom_data, 'aeiu', aeiu )
473            CALL iom_get  ( inum, jpdom_data, 'aeiv', aeiv )
474            CALL iom_close( inum )
475            !
476         CASE(  30  )       !==  fixed 3D shape  ==!
477            IF(lwp) WRITE(numout,*) '          tracer mixing coef. = F( latitude, longitude, depth )'
478            CALL ldf_c2d( 'TRA', 'LAP', rn_aeiv_0, aeiu, aeiv )    ! surface value proportional to scale factor
479            !                                                 ! reduction with depth
480            CALL ldf_c1d( 'TRA', r1_4, aeiu(:,:,1), aeiv(:,:,1), aeiu, aeiv )
481            !
482         CASE DEFAULT
483            CALL ctl_stop('ldf_tra_init: wrong choice for nn_aei_ijk_t, the type of space-time variation of aei')
484         END SELECT
485         !
486      ELSE
487          IF(lwp) WRITE(numout,*) '   eddy induced velocity param is NOT used neither diagnosed'
488          ln_ldfeiv_dia = .FALSE.
489      ENDIF
490      !                   
491   END SUBROUTINE ldf_eiv_init
492
493
494   SUBROUTINE ldf_eiv( kt, paei0, paeiu, paeiv )
495      !!----------------------------------------------------------------------
496      !!                  ***  ROUTINE ldf_eiv  ***
497      !!
498      !! ** Purpose :   Compute the eddy induced velocity coefficient from the
499      !!              growth rate of baroclinic instability.
500      !!
501      !! ** Method  :   coefficient function of the growth rate of baroclinic instability
502      !!
503      !! Reference : Treguier et al. JPO 1997   ; Held and Larichev JAS 1996
504      !!----------------------------------------------------------------------
505      INTEGER                         , INTENT(in   ) ::   kt             ! ocean time-step index
506      REAL(wp)                        , INTENT(inout) ::   paei0          ! max value            [m2/s]
507      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   paeiu, paeiv   ! eiv coefficient      [m2/s]
508      !
509      INTEGER  ::   ji, jj, jk    ! dummy loop indices
510      REAL(wp) ::   zfw, ze3w, zn2, z1_f20, zaht, zaht_min, zzaei   ! local scalars
511      REAL(wp), DIMENSION(:,:), POINTER ::   zn, zah, zhw, zross, zaeiw   ! 2D workspace
512      !!----------------------------------------------------------------------
513      !
514      IF( nn_timing == 1 )   CALL timing_start('ldf_eiv')
515      !
516      CALL wrk_alloc( jpi,jpj,   zn, zah, zhw, zross, zaeiw )
517      !     
518      zn   (:,:) = 0._wp      ! Local initialization
519      zhw  (:,:) = 5._wp
520      zah  (:,:) = 0._wp
521      zross(:,:) = 0._wp
522      !                       ! Compute lateral diffusive coefficient at T-point
523      IF( ln_traldf_triad ) THEN
524         DO jk = 1, jpk
525            DO jj = 2, jpjm1
526               DO ji = 2, jpim1
527                  ! Take the max of N^2 and zero then take the vertical sum
528                  ! of the square root of the resulting N^2 ( required to compute
529                  ! internal Rossby radius Ro = .5 * sum_jpk(N) / f
530                  zn2 = MAX( rn2b(ji,jj,jk), 0._wp )
531                  zn(ji,jj) = zn(ji,jj) + SQRT( zn2 ) * e3w_n(ji,jj,jk)
532                  ! Compute elements required for the inverse time scale of baroclinic
533                  ! eddies using the isopycnal slopes calculated in ldfslp.F :
534                  ! T^-1 = sqrt(m_jpk(N^2*(r1^2+r2^2)*e3w))
535                  ze3w = e3w_n(ji,jj,jk) * tmask(ji,jj,jk)
536                  zah(ji,jj) = zah(ji,jj) + zn2 * wslp2(ji,jj,jk) * ze3w
537                  zhw(ji,jj) = zhw(ji,jj) + ze3w
538               END DO
539            END DO
540         END DO
541      ELSE
542         DO jk = 1, jpk
543            DO jj = 2, jpjm1
544               DO ji = 2, jpim1
545                  ! Take the max of N^2 and zero then take the vertical sum
546                  ! of the square root of the resulting N^2 ( required to compute
547                  ! internal Rossby radius Ro = .5 * sum_jpk(N) / f
548                  zn2 = MAX( rn2b(ji,jj,jk), 0._wp )
549                  zn(ji,jj) = zn(ji,jj) + SQRT( zn2 ) * e3w_n(ji,jj,jk)
550                  ! Compute elements required for the inverse time scale of baroclinic
551                  ! eddies using the isopycnal slopes calculated in ldfslp.F :
552                  ! T^-1 = sqrt(m_jpk(N^2*(r1^2+r2^2)*e3w))
553                  ze3w = e3w_n(ji,jj,jk) * tmask(ji,jj,jk)
554                  zah(ji,jj) = zah(ji,jj) + zn2 * ( wslpi(ji,jj,jk) * wslpi(ji,jj,jk)   &
555                     &                            + wslpj(ji,jj,jk) * wslpj(ji,jj,jk) ) * ze3w
556                  zhw(ji,jj) = zhw(ji,jj) + ze3w
557               END DO
558            END DO
559         END DO
560      END IF
561
562      DO jj = 2, jpjm1
563         DO ji = fs_2, fs_jpim1   ! vector opt.
564            zfw = MAX( ABS( 2. * omega * SIN( rad * gphit(ji,jj) ) ) , 1.e-10 )
565            ! Rossby radius at w-point taken < 40km and  > 2km
566            zross(ji,jj) = MAX( MIN( .4 * zn(ji,jj) / zfw, 40.e3 ), 2.e3 )
567            ! Compute aeiw by multiplying Ro^2 and T^-1
568            zaeiw(ji,jj) = zross(ji,jj) * zross(ji,jj) * SQRT( zah(ji,jj) / zhw(ji,jj) ) * tmask(ji,jj,1)
569         END DO
570      END DO
571
572!!gm      IF( cp_cfg == "orca" .AND. jp_cfg == 2 ) THEN   ! ORCA R2
573!!gm         DO jj = 2, jpjm1
574!!gm            DO ji = fs_2, fs_jpim1   ! vector opt.
575!!gm               ! Take the minimum between aeiw and 1000 m2/s over shelves (depth shallower than 650 m)
576!!gm               IF( mbkt(ji,jj) <= 20 )   zaeiw(ji,jj) = MIN( zaeiw(ji,jj), 1000. )
577!!gm            END DO
578!!gm         END DO
579!!gm      ENDIF
580
581      !                                         !==  Bound on eiv coeff.  ==!
582      z1_f20 = 1._wp / (  2._wp * omega * sin( rad * 20._wp )  )
583      DO jj = 2, jpjm1
584         DO ji = fs_2, fs_jpim1   ! vector opt.
585            zzaei = MIN( 1._wp, ABS( ff(ji,jj) * z1_f20 ) ) * zaeiw(ji,jj)       ! tropical decrease
586            zaeiw(ji,jj) = MIN( zzaei , paei0 )                                  ! Max value = paei0
587         END DO
588      END DO
589      CALL lbc_lnk( zaeiw(:,:), 'W', 1. )       ! lateral boundary condition
590      !               
591      DO jj = 2, jpjm1                          !== aei at u- and v-points  ==!
592         DO ji = fs_2, fs_jpim1   ! vector opt.
593            paeiu(ji,jj,1) = 0.5_wp * ( zaeiw(ji,jj) + zaeiw(ji+1,jj  ) ) * umask(ji,jj,1)
594            paeiv(ji,jj,1) = 0.5_wp * ( zaeiw(ji,jj) + zaeiw(ji  ,jj+1) ) * vmask(ji,jj,1)
595         END DO
596      END DO
597      CALL lbc_lnk( paeiu(:,:,1), 'U', 1. )   ;   CALL lbc_lnk( paeiv(:,:,1), 'V', 1. )      ! lateral boundary condition
598
599      DO jk = 2, jpkm1                          !==  deeper values equal the surface one  ==!
600         paeiu(:,:,jk) = paeiu(:,:,1) * umask(:,:,jk)
601         paeiv(:,:,jk) = paeiv(:,:,1) * vmask(:,:,jk)
602      END DO
603     
604      CALL wrk_dealloc( jpi,jpj,   zn, zah, zhw, zross, zaeiw )
605      !
606      IF( nn_timing == 1 )   CALL timing_stop('ldf_eiv')
607      !
608   END SUBROUTINE ldf_eiv
609
610
611   SUBROUTINE ldf_eiv_trp( kt, kit000, pun, pvn, pwn, cdtype )
612      !!----------------------------------------------------------------------
613      !!                  ***  ROUTINE ldf_eiv_trp  ***
614      !!
615      !! ** Purpose :   add to the input ocean transport the contribution of
616      !!              the eddy induced velocity parametrization.
617      !!
618      !! ** Method  :   The eddy induced transport is computed from a flux stream-
619      !!              function which depends on the slope of iso-neutral surfaces
620      !!              (see ldf_slp). For example, in the i-k plan :
621      !!                   psi_uw = mk(aeiu) e2u mi(wslpi)   [in m3/s]
622      !!                   Utr_eiv = - dk[psi_uw]
623      !!                   Vtr_eiv = + di[psi_uw]
624      !!                ln_ldfeiv_dia = T : output the associated streamfunction,
625      !!                                    velocity and heat transport (call ldf_eiv_dia)
626      !!
627      !! ** Action  : pun, pvn increased by the eiv transport
628      !!----------------------------------------------------------------------
629      INTEGER                         , INTENT(in   ) ::   kt       ! ocean time-step index
630      INTEGER                         , INTENT(in   ) ::   kit000   ! first time step index
631      CHARACTER(len=3)                , INTENT(in   ) ::   cdtype   ! =TRA or TRC (tracer indicator)
632      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pun      ! in : 3 ocean transport components   [m3/s]
633      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pvn      ! out: 3 ocean transport components   [m3/s]
634      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pwn      ! increased by the eiv                [m3/s]
635      !!
636      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
637      REAL(wp) ::   zuwk, zuwk1, zuwi, zuwi1   ! local scalars
638      REAL(wp) ::   zvwk, zvwk1, zvwj, zvwj1   !   -      -
639      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zpsi_uw, zpsi_vw
640      !!----------------------------------------------------------------------
641      !
642      IF( nn_timing == 1 )   CALL timing_start( 'ldf_eiv_trp')
643      !
644      CALL wrk_alloc( jpi,jpj,jpk,   zpsi_uw, zpsi_vw )
645
646      IF( kt == kit000 )  THEN
647         IF(lwp) WRITE(numout,*)
648         IF(lwp) WRITE(numout,*) 'ldf_eiv_trp : eddy induced advection on ', cdtype,' :'
649         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   add to velocity fields the eiv component'
650      ENDIF
651
652     
653      zpsi_uw(:,:, 1 ) = 0._wp   ;   zpsi_vw(:,:, 1 ) = 0._wp
654      zpsi_uw(:,:,jpk) = 0._wp   ;   zpsi_vw(:,:,jpk) = 0._wp
655      !
656      DO jk = 2, jpkm1
657         DO jj = 1, jpjm1
658            DO ji = 1, fs_jpim1   ! vector opt.
659               zpsi_uw(ji,jj,jk) = - 0.25_wp * e2u(ji,jj) * ( wslpi(ji,jj,jk  ) + wslpi(ji+1,jj,jk) )   &
660                  &                                       * ( aeiu (ji,jj,jk-1) + aeiu (ji  ,jj,jk) ) * umask(ji,jj,jk)
661               zpsi_vw(ji,jj,jk) = - 0.25_wp * e1v(ji,jj) * ( wslpj(ji,jj,jk  ) + wslpj(ji,jj+1,jk) )   &
662                  &                                       * ( aeiv (ji,jj,jk-1) + aeiv (ji,jj  ,jk) ) * vmask(ji,jj,jk)
663            END DO
664         END DO
665      END DO
666      !
667      DO jk = 1, jpkm1
668         DO jj = 1, jpjm1
669            DO ji = 1, fs_jpim1   ! vector opt.               
670               pun(ji,jj,jk) = pun(ji,jj,jk) - ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji,jj,jk+1) )
671               pvn(ji,jj,jk) = pvn(ji,jj,jk) - ( zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj,jk+1) )
672            END DO
673         END DO
674      END DO
675      DO jk = 1, jpkm1
676         DO jj = 2, jpjm1
677            DO ji = fs_2, fs_jpim1   ! vector opt.
678               pwn(ji,jj,jk) = pwn(ji,jj,jk) + (  zpsi_uw(ji,jj,jk) - zpsi_uw(ji-1,jj  ,jk)   &
679                  &                             + zpsi_vw(ji,jj,jk) - zpsi_vw(ji  ,jj-1,jk) )
680            END DO
681         END DO
682      END DO
683      !
684      !                              ! diagnose the eddy induced velocity and associated heat transport
685      IF( ln_ldfeiv_dia .AND. cdtype == 'TRA' )   CALL ldf_eiv_dia( zpsi_uw, zpsi_vw )
686      !
687      CALL wrk_dealloc( jpi,jpj,jpk,   zpsi_uw, zpsi_vw )
688      !
689      IF( nn_timing == 1 )   CALL timing_stop( 'ldf_eiv_trp')
690      !
691    END SUBROUTINE ldf_eiv_trp
692
693
694   SUBROUTINE ldf_eiv_dia( psi_uw, psi_vw )
695      !!----------------------------------------------------------------------
696      !!                  ***  ROUTINE ldf_eiv_dia  ***
697      !!
698      !! ** Purpose :   diagnose the eddy induced velocity and its associated
699      !!              vertically integrated heat transport.
700      !!
701      !! ** Method :
702      !!
703      !!----------------------------------------------------------------------
704      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   psi_uw, psi_vw   ! streamfunction   [m3/s]
705      !
706      INTEGER  ::   ji, jj, jk    ! dummy loop indices
707      REAL(wp) ::   zztmp   ! local scalar
708      REAL(wp), DIMENSION(:,:)  , POINTER ::   zw2d   ! 2D workspace
709      REAL(wp), DIMENSION(:,:,:), POINTER ::   zw3d   ! 3D workspace
710      !!----------------------------------------------------------------------
711      !
712      IF( nn_timing == 1 )  CALL timing_start( 'ldf_eiv_dia')
713      !
714      !                                                  !==  eiv stream function: output  ==!
715      CALL lbc_lnk( psi_uw, 'U', -1. )                         ! lateral boundary condition
716      CALL lbc_lnk( psi_vw, 'V', -1. )
717      !
718!!gm      CALL iom_put( "psi_eiv_uw", psi_uw )                 ! output
719!!gm      CALL iom_put( "psi_eiv_vw", psi_vw )
720      !
721      !                                                  !==  eiv velocities: calculate and output  ==!
722      CALL wrk_alloc( jpi,jpj,jpk,   zw3d )
723      !
724      zw3d(:,:,jpk) = 0._wp                                    ! bottom value always 0
725      !
726      DO jk = 1, jpkm1                                         ! e2u e3u u_eiv = -dk[psi_uw]
727         zw3d(:,:,jk) = ( psi_uw(:,:,jk+1) - psi_uw(:,:,jk) ) / ( e2u(:,:) * e3u_n(:,:,jk) )
728      END DO
729      CALL iom_put( "uoce_eiv", zw3d )
730      !
731      DO jk = 1, jpkm1                                         ! e1v e3v v_eiv = -dk[psi_vw]
732         zw3d(:,:,jk) = ( psi_vw(:,:,jk+1) - psi_vw(:,:,jk) ) / ( e1v(:,:) * e3v_n(:,:,jk) )
733      END DO
734      CALL iom_put( "voce_eiv", zw3d )
735      !
736      DO jk = 1, jpkm1                                         ! e1 e2 w_eiv = dk[psix] + dk[psix]
737         DO jj = 2, jpjm1
738            DO ji = fs_2, fs_jpim1  ! vector opt.
739               zw3d(ji,jj,jk) = (  psi_vw(ji,jj,jk) - psi_vw(ji  ,jj-1,jk)    &
740                  &              + psi_uw(ji,jj,jk) - psi_uw(ji-1,jj  ,jk)  ) / e1e2t(ji,jj)
741            END DO
742         END DO
743      END DO
744      CALL lbc_lnk( zw3d, 'T', 1. )      ! lateral boundary condition
745      CALL iom_put( "woce_eiv", zw3d )
746      !
747      CALL wrk_dealloc( jpi,jpj,jpk,   zw3d )
748      !     
749      !
750      IF( lk_diaar5 ) THEN                               !==  eiv heat transport: calculate and output  ==!
751         CALL wrk_alloc( jpi,jpj,   zw2d )
752         !
753         zztmp = 0.5_wp * rau0 * rcp 
754         zw2d(:,:) = 0._wp 
755         DO jk = 1, jpkm1
756            DO jj = 2, jpjm1
757               DO ji = fs_2, fs_jpim1   ! vector opt.
758                  zw2d(ji,jj) = zw2d(ji,jj) + zztmp * ( psi_uw(ji,jj,jk+1)      - psi_uw(ji,jj,jk)          )   &
759                     &                              * ( tsn   (ji,jj,jk,jp_tem) + tsn   (ji+1,jj,jk,jp_tem) ) 
760               END DO
761            END DO
762         END DO
763         CALL lbc_lnk( zw2d, 'U', -1. )
764         CALL iom_put( "ueiv_heattr", zw2d )                  ! heat transport in i-direction
765         zw2d(:,:) = 0._wp 
766         DO jk = 1, jpkm1
767            DO jj = 2, jpjm1
768               DO ji = fs_2, fs_jpim1   ! vector opt.
769                  zw2d(ji,jj) = zw2d(ji,jj) + zztmp * ( psi_vw(ji,jj,jk+1)      - psi_vw(ji,jj,jk)          )   &
770                     &                              * ( tsn   (ji,jj,jk,jp_tem) + tsn   (ji,jj+1,jk,jp_tem) ) 
771               END DO
772            END DO
773         END DO
774         CALL lbc_lnk( zw2d, 'V', -1. )
775         CALL iom_put( "veiv_heattr", zw2d )                  !  heat transport in i-direction
776         !
777         CALL wrk_dealloc( jpi,jpj,   zw2d )
778      ENDIF
779      !
780      IF( nn_timing == 1 )  CALL timing_stop( 'ldf_eiv_dia')     
781      !
782   END SUBROUTINE ldf_eiv_dia
783
784   !!======================================================================
785END MODULE ldftra
Note: See TracBrowser for help on using the repository browser.