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

source: trunk/NEMOGCM/NEMO/OPA_SRC/LDF/ldftra.F90 @ 6140

Last change on this file since 6140 was 6140, checked in by timgraham, 8 years ago

Merge of branches/2015/dev_merge_2015 back into trunk. Merge excludes NEMOGCM/TOOLS/OBSTOOLS/ for now due to issues with the change of file type. Will sort these manually with further commits.

Branch merged as follows:
In the working copy of branch ran:
svn merge svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk@HEAD
Small conflicts due to bug fixes applied to trunk since the dev_merge_2015 was copied. Bug fixes were applied to the branch as well so these were easy to resolve.
Branch committed at this stage

In working copy run:
svn switch svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk
to switch working copy

Run:
svn merge --reintegrate svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/branches/2015/dev_merge_2015
to merge the branch into the trunk and then commit - no conflicts at this stage.

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