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 @ 9124

Last change on this file since 9124 was 9124, checked in by gm, 6 years ago

dev_merge_2017: ln_timing instead of nn_timing + restricted timing to nemo_init and routine called by step in OPA_SRC

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