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/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/LDF – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/LDF/ldftra.F90 @ 9094

Last change on this file since 9094 was 9094, checked in by cetlod, 6 years ago

Use of lbclnk_multi in subdir LDF & TRA

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