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 NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/OCE/LDF – NEMO

source: NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/OCE/LDF/ldftra.F90 @ 13630

Last change on this file since 13630 was 13630, checked in by mocavero, 3 years ago

Add neighborhood collectives calls in the NEMO src - ticket #2496

  • Property svn:keywords set to Id
File size: 53.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 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_OFF       !: 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   !                                    != iso-neutral options =!
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   !                                            !  time invariant coefficients:  aht_0 = 1/2  Ud*Ld   (lap case)
62   !                                            !                                bht_0 = 1/12 Ud*Ld^3 (blp case)
63   REAL(wp), PUBLIC ::      rn_Ud               !: lateral diffusive velocity  [m/s]
64   REAL(wp), PUBLIC ::      rn_Ld               !: lateral diffusive length    [m]
65
66   !                                   !!* Namelist namtra_eiv : eddy induced velocity param. *
67   !                                    != Use/diagnose eiv =!
68   LOGICAL , PUBLIC ::   ln_ldfeiv           !: eddy induced velocity flag
69   LOGICAL , PUBLIC ::   ln_ldfeiv_dia       !: diagnose & output eiv streamfunction and velocity (IOM)
70   !                                    != Coefficients =!
71   INTEGER , PUBLIC ::   nn_aei_ijk_t        !: choice of time/space variation of the eiv coeff.
72   REAL(wp), PUBLIC ::      rn_Ue               !: lateral diffusive velocity  [m/s]
73   REAL(wp), PUBLIC ::      rn_Le               !: lateral diffusive length    [m]
74   
75   !                                  ! Flag to control the type of lateral diffusive operator
76   INTEGER, PARAMETER, PUBLIC ::   np_ERROR  =-10   ! error in specification of lateral diffusion
77   INTEGER, PARAMETER, PUBLIC ::   np_no_ldf = 00   ! without operator (i.e. no lateral diffusive trend)
78   !                          !!      laplacian     !    bilaplacian    !
79   INTEGER, PARAMETER, PUBLIC ::   np_lap    = 10   ,   np_blp    = 20  ! iso-level operator
80   INTEGER, PARAMETER, PUBLIC ::   np_lap_i  = 11   ,   np_blp_i  = 21  ! standard iso-neutral or geopotential operator
81   INTEGER, PARAMETER, PUBLIC ::   np_lap_it = 12   ,   np_blp_it = 22  ! triad    iso-neutral or geopotential operator
82
83   INTEGER , PUBLIC ::   nldf_tra      = 0         !: type of lateral diffusion used defined from ln_traldf_... (namlist logicals)
84   LOGICAL , PUBLIC ::   l_ldftra_time = .FALSE.   !: flag for time variation of the lateral eddy diffusivity coef.
85   LOGICAL , PUBLIC ::   l_ldfeiv_time = .FALSE.   !: flag for time variation of the eiv coef.
86
87   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ahtu, ahtv   !: eddy diffusivity coef. at U- and V-points   [m2/s]
88   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   aeiu, aeiv   !: eddy induced velocity coeff.                [m2/s]
89
90   REAL(wp) ::   aht0, aei0               ! constant eddy coefficients (deduced from namelist values)                     [m2/s]
91   REAL(wp) ::   r1_2  = 0.5_wp           ! =1/2
92   REAL(wp) ::   r1_4  = 0.25_wp          ! =1/4
93   REAL(wp) ::   r1_12 = 1._wp / 12._wp   ! =1/12
94
95   !! * Substitutions
96#  include "do_loop_substitute.h90"
97#  include "domzgr_substitute.h90"
98   !!----------------------------------------------------------------------
99   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
100   !! $Id$
101   !! Software governed by the CeCILL license (see ./LICENSE)
102   !!----------------------------------------------------------------------
103CONTAINS
104
105   SUBROUTINE ldf_tra_init
106      !!----------------------------------------------------------------------
107      !!                  ***  ROUTINE ldf_tra_init  ***
108      !!
109      !! ** Purpose :   initializations of the tracer lateral mixing coeff.
110      !!
111      !! ** Method  : * the eddy diffusivity coef. specification depends on:
112      !!
113      !!    ln_traldf_lap = T     laplacian operator
114      !!    ln_traldf_blp = T   bilaplacian operator
115      !!
116      !!    nn_aht_ijk_t  =  0 => = constant
117      !!                  !
118      !!                  = 10 => = F(z) : constant with a reduction of 1/4 with depth
119      !!                  !
120      !!                  =-20 => = F(i,j)   = shape read in 'eddy_diffusivity.nc' file
121      !!                  = 20    = F(i,j)   = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
122      !!                  = 21    = F(i,j,t) = F(growth rate of baroclinic instability)
123      !!                  !
124      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_diffusivity.nc' file
125      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
126      !!                  = 31    = F(i,j,k,t) = F(local velocity) (  1/2  |u|e     laplacian operator
127      !!                                                           or 1/12 |u|e^3 bilaplacian operator )
128      !!              * initialisation of the eddy induced velocity coefficient by a call to ldf_eiv_init
129      !!           
130      !! ** action  : ahtu, ahtv initialized one for all or l_ldftra_time set to true
131      !!              aeiu, aeiv initialized one for all or l_ldfeiv_time set to true
132      !!----------------------------------------------------------------------
133      INTEGER  ::   jk                             ! dummy loop indices
134      INTEGER  ::   ioptio, ierr, inum, ios, inn   ! local integer
135      REAL(wp) ::   zah_max, zUfac                 !   -      -
136      CHARACTER(len=5) ::   cl_Units               ! units (m2/s or m4/s)
137      !!
138      NAMELIST/namtra_ldf/ ln_traldf_OFF, ln_traldf_lap  , ln_traldf_blp  ,   &   ! type of operator
139         &                 ln_traldf_lev, ln_traldf_hor  , ln_traldf_triad,   &   ! acting direction of the operator
140         &                 ln_traldf_iso, ln_traldf_msc  ,  rn_slpmax     ,   &   ! option for iso-neutral operator
141         &                 ln_triad_iso , ln_botmix_triad, rn_sw_triad    ,   &   ! option for triad operator
142         &                 nn_aht_ijk_t , rn_Ud          , rn_Ld                  ! lateral eddy coefficient
143      !!----------------------------------------------------------------------
144      !
145      IF(lwp) THEN                      ! control print
146         WRITE(numout,*)
147         WRITE(numout,*) 'ldf_tra_init : lateral tracer diffusion'
148         WRITE(numout,*) '~~~~~~~~~~~~ '
149      ENDIF
150     
151      !
152      !  Choice of lateral tracer physics
153      ! =================================
154      !
155      READ  ( numnam_ref, namtra_ldf, IOSTAT = ios, ERR = 901)
156901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtra_ldf in reference namelist' )
157      READ  ( numnam_cfg, namtra_ldf, IOSTAT = ios, ERR = 902 )
158902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtra_ldf in configuration namelist' )
159      IF(lwm) WRITE( numond, namtra_ldf )
160      !
161      IF(lwp) THEN                      ! control print
162         WRITE(numout,*) '   Namelist : namtra_ldf --- lateral mixing parameters (type, direction, coefficients)'
163         WRITE(numout,*) '      type :'
164         WRITE(numout,*) '         no explicit diffusion                   ln_traldf_OFF   = ', ln_traldf_OFF
165         WRITE(numout,*) '         laplacian operator                      ln_traldf_lap   = ', ln_traldf_lap
166         WRITE(numout,*) '         bilaplacian operator                    ln_traldf_blp   = ', ln_traldf_blp
167         WRITE(numout,*) '      direction of action :'
168         WRITE(numout,*) '         iso-level                               ln_traldf_lev   = ', ln_traldf_lev
169         WRITE(numout,*) '         horizontal (geopotential)               ln_traldf_hor   = ', ln_traldf_hor
170         WRITE(numout,*) '         iso-neutral Madec operator              ln_traldf_iso   = ', ln_traldf_iso
171         WRITE(numout,*) '         iso-neutral triad operator              ln_traldf_triad = ', ln_traldf_triad
172         WRITE(numout,*) '            use the Method of Stab. Correction   ln_traldf_msc   = ', ln_traldf_msc
173         WRITE(numout,*) '            maximum isoppycnal slope             rn_slpmax       = ', rn_slpmax
174         WRITE(numout,*) '            pure lateral mixing in ML            ln_triad_iso    = ', ln_triad_iso
175         WRITE(numout,*) '            switching triad or not               rn_sw_triad     = ', rn_sw_triad
176         WRITE(numout,*) '            lateral mixing on bottom             ln_botmix_triad = ', ln_botmix_triad
177         WRITE(numout,*) '      coefficients :'
178         WRITE(numout,*) '         type of time-space variation            nn_aht_ijk_t    = ', nn_aht_ijk_t
179         WRITE(numout,*) '            lateral diffusive velocity (if cst)  rn_Ud           = ', rn_Ud, ' m/s'
180         WRITE(numout,*) '            lateral diffusive length   (if cst)  rn_Ld           = ', rn_Ld, ' m'
181      ENDIF
182      !
183      !
184      ! Operator and its acting direction   (set nldf_tra) 
185      ! =================================
186      !
187      nldf_tra = np_ERROR
188      ioptio   = 0
189      IF( ln_traldf_OFF ) THEN   ;   nldf_tra = np_no_ldf   ;   ioptio = ioptio + 1   ;   ENDIF
190      IF( ln_traldf_lap ) THEN   ;                              ioptio = ioptio + 1   ;   ENDIF
191      IF( ln_traldf_blp ) THEN   ;                              ioptio = ioptio + 1   ;   ENDIF
192      IF( ioptio /=  1  )   CALL ctl_stop( 'tra_ldf_init: use ONE of the 3 operator options (NONE/lap/blp)' )
193      !
194      IF( .NOT.ln_traldf_OFF ) THEN    !==  direction ==>> type of operator  ==!
195         ioptio = 0
196         IF( ln_traldf_lev   )   ioptio = ioptio + 1
197         IF( ln_traldf_hor   )   ioptio = ioptio + 1
198         IF( ln_traldf_iso   )   ioptio = ioptio + 1
199         IF( ln_traldf_triad )   ioptio = ioptio + 1
200         IF( ioptio /=  1  )   CALL ctl_stop( 'tra_ldf_init: use ONE direction (level/hor/iso/triad)' )
201         !
202         !                                ! defined the type of lateral diffusion from ln_traldf_... logicals
203         ierr = 0
204         IF ( ln_traldf_lap ) THEN        ! laplacian operator
205            IF ( ln_zco ) THEN                  ! z-coordinate
206               IF ( ln_traldf_lev   )   nldf_tra = np_lap     ! iso-level = horizontal (no rotation)
207               IF ( ln_traldf_hor   )   nldf_tra = np_lap     ! iso-level = horizontal (no rotation)
208               IF ( ln_traldf_iso   )   nldf_tra = np_lap_i   ! iso-neutral: standard  (   rotation)
209               IF ( ln_traldf_triad )   nldf_tra = np_lap_it  ! iso-neutral: triad     (   rotation)
210            ENDIF
211            IF ( ln_zps ) THEN                  ! z-coordinate with partial step
212               IF ( ln_traldf_lev   )   ierr     = 1          ! iso-level not allowed
213               IF ( ln_traldf_hor   )   nldf_tra = np_lap     ! horizontal             (no rotation)
214               IF ( ln_traldf_iso   )   nldf_tra = np_lap_i   ! iso-neutral: standard     (rotation)
215               IF ( ln_traldf_triad )   nldf_tra = np_lap_it  ! iso-neutral: triad        (rotation)
216            ENDIF
217            IF ( ln_sco ) THEN                  ! s-coordinate
218               IF ( ln_traldf_lev   )   nldf_tra = np_lap     ! iso-level              (no rotation)
219               IF ( ln_traldf_hor   )   nldf_tra = np_lap_i   ! horizontal             (   rotation)
220               IF ( ln_traldf_iso   )   nldf_tra = np_lap_i   ! iso-neutral: standard  (   rotation)
221               IF ( ln_traldf_triad )   nldf_tra = np_lap_it  ! iso-neutral: triad     (   rotation)
222            ENDIF
223         ENDIF
224         !
225         IF( ln_traldf_blp ) THEN         ! bilaplacian operator
226            IF ( ln_zco ) THEN                  ! z-coordinate
227               IF ( ln_traldf_lev   )   nldf_tra = np_blp     ! iso-level = horizontal (no rotation)
228               IF ( ln_traldf_hor   )   nldf_tra = np_blp     ! iso-level = horizontal (no rotation)
229               IF ( ln_traldf_iso   )   nldf_tra = np_blp_i   ! iso-neutral: standard  (   rotation)
230               IF ( ln_traldf_triad )   nldf_tra = np_blp_it  ! iso-neutral: triad     (   rotation)
231            ENDIF
232            IF ( ln_zps ) THEN                  ! z-coordinate with partial step
233               IF ( ln_traldf_lev   )   ierr     = 1          ! iso-level not allowed
234               IF ( ln_traldf_hor   )   nldf_tra = np_blp     ! horizontal             (no rotation)
235               IF ( ln_traldf_iso   )   nldf_tra = np_blp_i   ! iso-neutral: standard  (   rotation)
236               IF ( ln_traldf_triad )   nldf_tra = np_blp_it  ! iso-neutral: triad     (   rotation)
237            ENDIF
238            IF ( ln_sco ) THEN                  ! s-coordinate
239               IF ( ln_traldf_lev   )   nldf_tra = np_blp     ! iso-level              (no rotation)
240               IF ( ln_traldf_hor   )   nldf_tra = np_blp_it  ! horizontal             (   rotation)
241               IF ( ln_traldf_iso   )   nldf_tra = np_blp_i   ! iso-neutral: standard  (   rotation)
242               IF ( ln_traldf_triad )   nldf_tra = np_blp_it  ! iso-neutral: triad     (   rotation)
243            ENDIF
244         ENDIF
245         IF ( ierr == 1 )   CALL ctl_stop( 'iso-level in z-partial step, not allowed' )
246      ENDIF
247      !
248      IF( ln_isfcav .AND. ln_traldf_triad )   CALL ctl_stop( ' ice shelf cavity and traldf_triad not tested' )
249           !
250      IF(  nldf_tra == np_lap_i .OR. nldf_tra == np_lap_it .OR. &
251         & nldf_tra == np_blp_i .OR. nldf_tra == np_blp_it  )   l_ldfslp = .TRUE.    ! slope of neutral surfaces required
252      !
253      IF( ln_traldf_blp .AND. ( ln_traldf_iso .OR. ln_traldf_triad) ) THEN     ! iso-neutral bilaplacian need MSC
254         IF( .NOT.ln_traldf_msc )   CALL ctl_stop( 'tra_ldf_init: iso-neutral bilaplacian requires ln_traldf_msc=.true.' )
255      ENDIF
256      !
257      IF(lwp) THEN
258         WRITE(numout,*)
259         SELECT CASE( nldf_tra )
260         CASE( np_no_ldf )   ;   WRITE(numout,*) '   ==>>>   NO lateral diffusion'
261         CASE( np_lap    )   ;   WRITE(numout,*) '   ==>>>   laplacian iso-level operator'
262         CASE( np_lap_i  )   ;   WRITE(numout,*) '   ==>>>   Rotated laplacian operator (standard)'
263         CASE( np_lap_it )   ;   WRITE(numout,*) '   ==>>>   Rotated laplacian operator (triad)'
264         CASE( np_blp    )   ;   WRITE(numout,*) '   ==>>>   bilaplacian iso-level operator'
265         CASE( np_blp_i  )   ;   WRITE(numout,*) '   ==>>>   Rotated bilaplacian operator (standard)'
266         CASE( np_blp_it )   ;   WRITE(numout,*) '   ==>>>   Rotated bilaplacian operator (triad)'
267         END SELECT
268         WRITE(numout,*)
269      ENDIF
270
271      !
272      !  Space/time variation of eddy coefficients
273      ! ===========================================
274      !
275      l_ldftra_time = .FALSE.                ! no time variation except in case defined below
276      !
277      IF( ln_traldf_OFF ) THEN               !== no explicit diffusive operator  ==!
278         !
279         IF(lwp) WRITE(numout,*) '   ==>>>   No diffusive operator selected. ahtu and ahtv are not allocated'
280         RETURN
281         !
282      ELSE                                   !==  a lateral diffusion operator is used  ==!
283         !
284         !                                         ! allocate the aht arrays
285         ALLOCATE( ahtu(jpi,jpj,jpk) , ahtv(jpi,jpj,jpk) , STAT=ierr )
286         IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_tra_init: failed to allocate arrays')
287         !
288         ahtu(:,:,jpk) = 0._wp                     ! last level always 0 
289         ahtv(:,:,jpk) = 0._wp
290         !.
291         !                                         ! value of lap/blp eddy mixing coef.
292         IF(     ln_traldf_lap ) THEN   ;   zUfac = r1_2 *rn_Ud   ;   inn = 1   ;   cl_Units = ' m2/s'   !   laplacian
293         ELSEIF( ln_traldf_blp ) THEN   ;   zUfac = r1_12*rn_Ud   ;   inn = 3   ;   cl_Units = ' m4/s'   ! bilaplacian
294         ENDIF
295         aht0    = zUfac *    rn_Ld**inn              ! mixing coefficient
296         zah_max = zUfac * (ra*rad)**inn              ! maximum reachable coefficient (value at the Equator for e1=1 degree)
297         !
298         !
299         SELECT CASE(  nn_aht_ijk_t  )             !* Specification of space-time variations of ahtu, ahtv
300         !
301         CASE(   0  )      !==  constant  ==!
302            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = constant = ', aht0, cl_Units
303            ahtu(:,:,1:jpkm1) = aht0
304            ahtv(:,:,1:jpkm1) = aht0
305            !
306         CASE(  10  )      !==  fixed profile  ==!
307            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( depth )'
308            IF(lwp) WRITE(numout,*) '           surface eddy diffusivity = constant = ', aht0, cl_Units
309            ahtu(:,:,1) = aht0                        ! constant surface value
310            ahtv(:,:,1) = aht0
311            CALL ldf_c1d( 'TRA', ahtu(:,:,1), ahtv(:,:,1), ahtu, ahtv )
312            !
313         CASE ( -20 )      !== fixed horizontal shape and magnitude read in file  ==!
314            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F(i,j) read in eddy_diffusivity.nc file'
315            CALL iom_open( 'eddy_diffusivity_2D.nc', inum )
316            CALL iom_get ( inum, jpdom_global, 'ahtu_2D', ahtu(:,:,1), cd_type = 'U', psgn = 1._wp )
317            CALL iom_get ( inum, jpdom_global, 'ahtv_2D', ahtv(:,:,1), cd_type = 'V', psgn = 1._wp )
318            CALL iom_close( inum )
319            DO jk = 2, jpkm1
320               ahtu(:,:,jk) = ahtu(:,:,1)
321               ahtv(:,:,jk) = ahtv(:,:,1)
322            END DO
323            !
324         CASE(  20  )      !== fixed horizontal shape  ==!
325            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( e1, e2 ) or F( e1^3, e2^3 ) (lap or blp case)'
326            IF(lwp) WRITE(numout,*) '           using a fixed diffusive velocity = ', rn_Ud,' m/s   and   Ld = Max(e1,e2)'
327            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
328            CALL ldf_c2d( 'TRA', zUfac      , inn        , ahtu, ahtv )    ! value proportional to scale factor^inn
329            !
330         CASE(  21  )      !==  time varying 2D field  ==!
331            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( latitude, longitude, time )'
332            IF(lwp) WRITE(numout,*) '                            = F( growth rate of baroclinic instability )'
333            IF(lwp) WRITE(numout,*) '            min value = 0.2 * aht0 (with aht0= 1/2 rn_Ud*rn_Ld)'
334            IF(lwp) WRITE(numout,*) '            max value =       aei0 (with aei0=1/2 rn_Ue*Le  increased to aht0 within 20N-20S'
335            !
336            l_ldftra_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
337            !
338            IF( ln_traldf_blp )   CALL ctl_stop( 'ldf_tra_init: aht=F( growth rate of baroc. insta .)',   &
339               &                                 '              incompatible with bilaplacian operator' )
340            !
341         CASE( -30  )      !== fixed 3D shape read in file  ==!
342            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F(i,j,k) read in eddy_diffusivity.nc file'
343            CALL iom_open( 'eddy_diffusivity_3D.nc', inum )
344            CALL iom_get ( inum, jpdom_global, 'ahtu_3D', ahtu, cd_type = 'U', psgn = 1._wp )
345            CALL iom_get ( inum, jpdom_global, 'ahtv_3D', ahtv, cd_type = 'V', psgn = 1._wp )
346            CALL iom_close( inum )
347            !
348         CASE(  30  )      !==  fixed 3D shape  ==!
349            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( latitude, longitude, depth )'
350            IF(lwp) WRITE(numout,*) '           using a fixed diffusive velocity = ', rn_Ud,' m/s   and   Ld = Max(e1,e2)'
351            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
352            CALL ldf_c2d( 'TRA', zUfac      , inn        , ahtu, ahtv )    ! surface value proportional to scale factor^inn
353            CALL ldf_c1d( 'TRA', ahtu(:,:,1), ahtv(:,:,1), ahtu, ahtv )    ! reduction with depth
354            !
355         CASE(  31  )      !==  time varying 3D field  ==!
356            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( latitude, longitude, depth , time )'
357            IF(lwp) WRITE(numout,*) '           proportional to the velocity : 1/2 |u|e or 1/12 |u|e^3'
358            !
359            l_ldftra_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
360            !
361         CASE DEFAULT
362            CALL ctl_stop('ldf_tra_init: wrong choice for nn_aht_ijk_t, the type of space-time variation of aht')
363         END SELECT
364         !
365         IF( .NOT.l_ldftra_time ) THEN             !* No time variation
366            IF(     ln_traldf_lap ) THEN                 !   laplacian operator (mask only)
367               ahtu(:,:,1:jpkm1) =       ahtu(:,:,1:jpkm1)   * umask(:,:,1:jpkm1)
368               ahtv(:,:,1:jpkm1) =       ahtv(:,:,1:jpkm1)   * vmask(:,:,1:jpkm1)
369            ELSEIF( ln_traldf_blp ) THEN                 ! bilaplacian operator (square root + mask)
370               ahtu(:,:,1:jpkm1) = SQRT( ahtu(:,:,1:jpkm1) ) * umask(:,:,1:jpkm1)
371               ahtv(:,:,1:jpkm1) = SQRT( ahtv(:,:,1:jpkm1) ) * vmask(:,:,1:jpkm1)
372            ENDIF
373         ENDIF
374         !
375      ENDIF
376      !
377   END SUBROUTINE ldf_tra_init
378
379
380   SUBROUTINE ldf_tra( kt, Kbb, Kmm )
381      !!----------------------------------------------------------------------
382      !!                  ***  ROUTINE ldf_tra  ***
383      !!
384      !! ** Purpose :   update at kt the tracer lateral mixing coeff. (aht and aeiv)
385      !!
386      !! ** Method  : * time varying eddy diffusivity coefficients:
387      !!
388      !!    nn_aei_ijk_t = 21    aeiu, aeiv = F(i,j,  t) = F(growth rate of baroclinic instability)
389      !!                                                   with a reduction to 0 in vicinity of the Equator
390      !!    nn_aht_ijk_t = 21    ahtu, ahtv = F(i,j,  t) = F(growth rate of baroclinic instability)
391      !!
392      !!                 = 31    ahtu, ahtv = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
393      !!                                                                     or |u|e^3/12 bilaplacian operator )
394      !!
395      !!              * time varying EIV coefficients: call to ldf_eiv routine
396      !!
397      !! ** action  :   ahtu, ahtv   update at each time step   
398      !!                aeiu, aeiv      -       -     -    -   (if ln_ldfeiv=T)
399      !!----------------------------------------------------------------------
400      INTEGER, INTENT(in) ::   kt   ! time step
401      INTEGER, INTENT(in) ::   Kbb, Kmm   ! ocean time level indices
402      !
403      INTEGER  ::   ji, jj, jk   ! dummy loop indices
404      REAL(wp) ::   zaht, zahf, zaht_min, zDaht, z1_f20   ! local scalar
405      !!----------------------------------------------------------------------
406      !
407      IF( ln_ldfeiv .AND. nn_aei_ijk_t == 21 ) THEN       ! eddy induced velocity coefficients
408         !                                ! =F(growth rate of baroclinic instability)
409         !                                ! max value aeiv_0 ; decreased to 0 within 20N-20S
410         CALL ldf_eiv( kt, aei0, aeiu, aeiv, Kmm )
411      ENDIF
412      !
413      SELECT CASE(  nn_aht_ijk_t  )       ! Eddy diffusivity coefficients
414      !
415      CASE(  21  )       !==  time varying 2D field  ==!   = F( growth rate of baroclinic instability )
416         !                                             !   min value 0.2*aht0
417         !                                             !   max value aht0 (aei0 if nn_aei_ijk_t=21)
418         !                                             !   increase to aht0 within 20N-20S
419         IF( ln_ldfeiv .AND. nn_aei_ijk_t == 21 ) THEN   ! use the already computed aei.
420            ahtu(:,:,1) = aeiu(:,:,1)
421            ahtv(:,:,1) = aeiv(:,:,1)
422         ELSE                                            ! compute aht.
423            CALL ldf_eiv( kt, aht0, ahtu, ahtv, Kmm )
424         ENDIF
425         !
426         z1_f20   = 1._wp / (  2._wp * omega * SIN( rad * 20._wp )  )   ! 1 / ff(20 degrees)   
427         zaht_min = 0.2_wp * aht0                                       ! minimum value for aht
428         zDaht    = aht0 - zaht_min                                     
429         DO_2D( 1, 1, 1, 1 )
430            !!gm CAUTION : here we assume lat/lon grid in 20deg N/S band (like all ORCA cfg)
431            !!     ==>>>   The Coriolis value is identical for t- & u_points, and for v- and f-points
432            zaht = ( 1._wp -  MIN( 1._wp , ABS( ff_t(ji,jj) * z1_f20 ) ) ) * zDaht
433            zahf = ( 1._wp -  MIN( 1._wp , ABS( ff_f(ji,jj) * z1_f20 ) ) ) * zDaht
434            ahtu(ji,jj,1) = (  MAX( zaht_min, ahtu(ji,jj,1) ) + zaht  )     ! min value zaht_min
435            ahtv(ji,jj,1) = (  MAX( zaht_min, ahtv(ji,jj,1) ) + zahf  )     ! increase within 20S-20N
436         END_2D
437         DO jk = 1, jpkm1                             ! deeper value = surface value + mask for all levels
438            ahtu(:,:,jk) = ahtu(:,:,1) * umask(:,:,jk)
439            ahtv(:,:,jk) = ahtv(:,:,1) * vmask(:,:,jk)
440         END DO
441         !
442      CASE(  31  )       !==  time varying 3D field  ==!   = F( local velocity )
443         IF( ln_traldf_lap     ) THEN          !   laplacian operator |u| e /12
444            DO jk = 1, jpkm1
445               ahtu(:,:,jk) = ABS( uu(:,:,jk,Kbb) ) * e1u(:,:) * r1_12   ! n.b. uu,vv are masked
446               ahtv(:,:,jk) = ABS( vv(:,:,jk,Kbb) ) * e2v(:,:) * r1_12
447            END DO
448         ELSEIF( ln_traldf_blp ) THEN      ! bilaplacian operator      sqrt( |u| e^3 /12 ) = sqrt( |u| e /12 ) * e
449            DO jk = 1, jpkm1
450               ahtu(:,:,jk) = SQRT(  ABS( uu(:,:,jk,Kbb) ) * e1u(:,:) * r1_12  ) * e1u(:,:)
451               ahtv(:,:,jk) = SQRT(  ABS( vv(:,:,jk,Kbb) ) * e2v(:,:) * r1_12  ) * e2v(:,:)
452            END DO
453         ENDIF
454         !
455      END SELECT
456      !
457      CALL iom_put( "ahtu_2d", ahtu(:,:,1) )   ! surface u-eddy diffusivity coeff.
458      CALL iom_put( "ahtv_2d", ahtv(:,:,1) )   ! surface v-eddy diffusivity coeff.
459      CALL iom_put( "ahtu_3d", ahtu(:,:,:) )   ! 3D      u-eddy diffusivity coeff.
460      CALL iom_put( "ahtv_3d", ahtv(:,:,:) )   ! 3D      v-eddy diffusivity coeff.
461      !
462      IF( ln_ldfeiv ) THEN
463        CALL iom_put( "aeiu_2d", aeiu(:,:,1) )   ! surface u-EIV coeff.
464        CALL iom_put( "aeiv_2d", aeiv(:,:,1) )   ! surface v-EIV coeff.
465        CALL iom_put( "aeiu_3d", aeiu(:,:,:) )   ! 3D      u-EIV coeff.
466        CALL iom_put( "aeiv_3d", aeiv(:,:,:) )   ! 3D      v-EIV coeff.
467      ENDIF
468      !
469   END SUBROUTINE ldf_tra
470
471
472   SUBROUTINE ldf_eiv_init
473      !!----------------------------------------------------------------------
474      !!                  ***  ROUTINE ldf_eiv_init  ***
475      !!
476      !! ** Purpose :   initialization of the eiv coeff. from namelist choices.
477      !!
478      !! ** Method  :   the eiv diffusivity coef. specification depends on:
479      !!    nn_aei_ijk_t  =  0 => = constant
480      !!                  !
481      !!                  = 10 => = F(z) : constant with a reduction of 1/4 with depth
482      !!                  !
483      !!                  =-20 => = F(i,j)   = shape read in 'eddy_diffusivity.nc' file
484      !!                  = 20    = F(i,j)   = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
485      !!                  = 21    = F(i,j,t) = F(growth rate of baroclinic instability)
486      !!                  !
487      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_diffusivity.nc' file
488      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
489      !!
490      !! ** Action  :   aeiu , aeiv   :  initialized one for all or l_ldftra_time set to true
491      !!                l_ldfeiv_time : =T if EIV coefficients vary with time
492      !!----------------------------------------------------------------------
493      INTEGER  ::   jk                     ! dummy loop indices
494      INTEGER  ::   ierr, inum, ios, inn   ! local integer
495      REAL(wp) ::   zah_max, zUfac         !   -   scalar
496      !!
497      NAMELIST/namtra_eiv/ ln_ldfeiv   , ln_ldfeiv_dia,   &   ! eddy induced velocity (eiv)
498         &                 nn_aei_ijk_t, rn_Ue, rn_Le         ! eiv  coefficient
499      !!----------------------------------------------------------------------
500      !
501      IF(lwp) THEN                      ! control print
502         WRITE(numout,*)
503         WRITE(numout,*) 'ldf_eiv_init : eddy induced velocity parametrization'
504         WRITE(numout,*) '~~~~~~~~~~~~ '
505      ENDIF
506      !
507      READ  ( numnam_ref, namtra_eiv, IOSTAT = ios, ERR = 901)
508901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtra_eiv in reference namelist' )
509      !
510      READ  ( numnam_cfg, namtra_eiv, IOSTAT = ios, ERR = 902 )
511902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtra_eiv in configuration namelist' )
512      IF(lwm)  WRITE ( numond, namtra_eiv )
513
514      IF(lwp) THEN                      ! control print
515         WRITE(numout,*) '   Namelist namtra_eiv : '
516         WRITE(numout,*) '      Eddy Induced Velocity (eiv) param.         ln_ldfeiv     = ', ln_ldfeiv
517         WRITE(numout,*) '      eiv streamfunction & velocity diag.        ln_ldfeiv_dia = ', ln_ldfeiv_dia
518         WRITE(numout,*) '      coefficients :'
519         WRITE(numout,*) '         type of time-space variation            nn_aei_ijk_t  = ', nn_aht_ijk_t
520         WRITE(numout,*) '         lateral diffusive velocity (if cst)     rn_Ue         = ', rn_Ue, ' m/s'
521         WRITE(numout,*) '         lateral diffusive length   (if cst)     rn_Le         = ', rn_Le, ' m'
522         WRITE(numout,*)
523      ENDIF
524      !
525      l_ldfeiv_time = .FALSE.       ! no time variation except in case defined below
526      !
527      !
528      IF( .NOT.ln_ldfeiv ) THEN     !== Parametrization not used  ==!
529         !
530         IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity param is NOT used'
531         ln_ldfeiv_dia = .FALSE.
532         !
533      ELSE                          !== use the parametrization  ==!
534         !
535         IF(lwp) WRITE(numout,*) '   ==>>>   use eddy induced velocity parametrization'
536         IF(lwp) WRITE(numout,*)
537         !
538         IF( ln_traldf_blp )   CALL ctl_stop( 'ldf_eiv_init: eddy induced velocity ONLY with laplacian diffusivity' )
539         !
540         IF( .NOT.( ln_traldf_iso .OR. ln_traldf_triad ) )   &
541           &                  CALL ctl_stop( 'ln_ldfeiv=T requires iso-neutral laplacian diffusion' )
542         !                                != allocate the aei arrays
543         ALLOCATE( aeiu(jpi,jpj,jpk), aeiv(jpi,jpj,jpk), STAT=ierr )
544         IF( ierr /= 0 )   CALL ctl_stop('STOP', 'ldf_eiv: failed to allocate arrays')
545         !
546         !                                != Specification of space-time variations of eaiu, aeiv
547         !
548         aeiu(:,:,jpk) = 0._wp               ! last level always 0 
549         aeiv(:,:,jpk) = 0._wp
550         !                                   ! value of EIV coef. (laplacian operator)
551         zUfac = r1_2 *rn_Ue                    ! velocity factor
552         inn = 1                                ! L-exponent
553         aei0    = zUfac *    rn_Le**inn        ! mixing coefficient
554         zah_max = zUfac * (ra*rad)**inn        ! maximum reachable coefficient (value at the Equator)
555
556         SELECT CASE( nn_aei_ijk_t )         !* Specification of space-time variations
557         !
558         CASE(   0  )                        !--  constant  --!
559            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = constant = ', aei0, ' m2/s'
560            aeiu(:,:,1:jpkm1) = aei0
561            aeiv(:,:,1:jpkm1) = aei0
562            !
563         CASE(  10  )                        !--  fixed profile  --!
564            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( depth )'
565            IF(lwp) WRITE(numout,*) '           surface eddy diffusivity = constant = ', aht0, ' m2/s'
566            aeiu(:,:,1) = aei0                  ! constant surface value
567            aeiv(:,:,1) = aei0
568            CALL ldf_c1d( 'TRA', aeiu(:,:,1), aeiv(:,:,1), aeiu, aeiv )
569            !
570         CASE ( -20 )                        !--  fixed horizontal shape read in file  --!
571            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F(i,j) read in eddy_diffusivity_2D.nc file'
572            CALL iom_open ( 'eddy_induced_velocity_2D.nc', inum )
573            CALL iom_get  ( inum, jpdom_global, 'aeiu', aeiu(:,:,1), cd_type = 'U', psgn = 1._wp )
574            CALL iom_get  ( inum, jpdom_global, 'aeiv', aeiv(:,:,1), cd_type = 'V', psgn = 1._wp )
575            CALL iom_close( inum )
576            DO jk = 2, jpkm1
577               aeiu(:,:,jk) = aeiu(:,:,1)
578               aeiv(:,:,jk) = aeiv(:,:,1)
579            END DO
580            !
581         CASE(  20  )                        !--  fixed horizontal shape  --!
582            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( e1, e2 )'
583            IF(lwp) WRITE(numout,*) '           using a fixed diffusive velocity = ', rn_Ue, ' m/s   and   Le = Max(e1,e2)'
584            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, ' m2/s   for e1=1°)'
585            CALL ldf_c2d( 'TRA', zUfac      , inn        , aeiu, aeiv )    ! value proportional to scale factor^inn
586            !
587         CASE(  21  )                        !--  time varying 2D field  --!
588            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( latitude, longitude, time )'
589            IF(lwp) WRITE(numout,*) '                                       = F( growth rate of baroclinic instability )'
590            IF(lwp) WRITE(numout,*) '           maximum allowed value: aei0 = ', aei0, ' m2/s'
591            !
592            l_ldfeiv_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
593            !
594         CASE( -30  )                        !-- fixed 3D shape read in file  --!
595            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F(i,j,k) read in eddy_diffusivity_3D.nc file'
596            CALL iom_open ( 'eddy_induced_velocity_3D.nc', inum )
597            CALL iom_get  ( inum, jpdom_global, 'aeiu', aeiu, cd_type = 'U', psgn = 1._wp )
598            CALL iom_get  ( inum, jpdom_global, 'aeiv', aeiv, cd_type = 'V', psgn = 1._wp )
599            CALL iom_close( inum )
600            !
601         CASE(  30  )                        !--  fixed 3D shape  --!
602            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( latitude, longitude, depth )'
603            CALL ldf_c2d( 'TRA', zUfac      , inn        , aeiu, aeiv )    ! surface value proportional to scale factor^inn
604            CALL ldf_c1d( 'TRA', aeiu(:,:,1), aeiv(:,:,1), aeiu, aeiv )    ! reduction with depth
605            !
606         CASE DEFAULT
607            CALL ctl_stop('ldf_tra_init: wrong choice for nn_aei_ijk_t, the type of space-time variation of aei')
608         END SELECT
609         !
610         IF( .NOT.l_ldfeiv_time ) THEN             !* mask if No time variation
611            DO jk = 1, jpkm1
612               aeiu(:,:,jk) = aeiu(:,:,jk) * umask(:,:,jk)
613               ahtv(:,:,jk) = ahtv(:,:,jk) * vmask(:,:,jk)
614            END DO
615         ENDIF
616         !
617      ENDIF
618      !                   
619   END SUBROUTINE ldf_eiv_init
620
621
622   SUBROUTINE ldf_eiv( kt, paei0, paeiu, paeiv, Kmm )
623      !!----------------------------------------------------------------------
624      !!                  ***  ROUTINE ldf_eiv  ***
625      !!
626      !! ** Purpose :   Compute the eddy induced velocity coefficient from the
627      !!              growth rate of baroclinic instability.
628      !!
629      !! ** Method  :   coefficient function of the growth rate of baroclinic instability
630      !!
631      !! Reference : Treguier et al. JPO 1997   ; Held and Larichev JAS 1996
632      !!----------------------------------------------------------------------
633      INTEGER                         , INTENT(in   ) ::   kt             ! ocean time-step index
634      INTEGER                         , INTENT(in   ) ::   Kmm            ! ocean time level indices
635      REAL(wp)                        , INTENT(inout) ::   paei0          ! max value            [m2/s]
636      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   paeiu, paeiv   ! eiv coefficient      [m2/s]
637      !
638      INTEGER  ::   ji, jj, jk    ! dummy loop indices
639      REAL(wp) ::   zfw, ze3w, zn2, z1_f20, zaht, zaht_min, zzaei    ! local scalars
640      REAL(wp), DIMENSION(jpi,jpj) ::   zn, zah, zhw, zRo, zaeiw   ! 2D workspace
641      !!----------------------------------------------------------------------
642      !
643      zn (:,:) = 0._wp        ! Local initialization
644      zhw(:,:) = 5._wp
645      zah(:,:) = 0._wp
646      zRo(:,:) = 0._wp
647      !                       ! Compute lateral diffusive coefficient at T-point
648      IF( ln_traldf_triad ) THEN
649         DO_3D( 0, 0, 0, 0, 1, jpk )
650            ! Take the max of N^2 and zero then take the vertical sum
651            ! of the square root of the resulting N^2 ( required to compute
652            ! internal Rossby radius Ro = .5 * sum_jpk(N) / f
653            zn2 = MAX( rn2b(ji,jj,jk), 0._wp )
654            zn(ji,jj) = zn(ji,jj) + SQRT( zn2 ) * e3w(ji,jj,jk,Kmm)
655            ! Compute elements required for the inverse time scale of baroclinic
656            ! eddies using the isopycnal slopes calculated in ldfslp.F :
657            ! T^-1 = sqrt(m_jpk(N^2*(r1^2+r2^2)*e3w))
658            ze3w = e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)
659            zah(ji,jj) = zah(ji,jj) + zn2 * wslp2(ji,jj,jk) * ze3w
660            zhw(ji,jj) = zhw(ji,jj) + ze3w
661         END_3D
662      ELSE
663         DO_3D( 0, 0, 0, 0, 1, jpk )
664            ! Take the max of N^2 and zero then take the vertical sum
665            ! of the square root of the resulting N^2 ( required to compute
666            ! internal Rossby radius Ro = .5 * sum_jpk(N) / f
667            zn2 = MAX( rn2b(ji,jj,jk), 0._wp )
668            zn(ji,jj) = zn(ji,jj) + SQRT( zn2 ) * e3w(ji,jj,jk,Kmm)
669            ! Compute elements required for the inverse time scale of baroclinic
670            ! eddies using the isopycnal slopes calculated in ldfslp.F :
671            ! T^-1 = sqrt(m_jpk(N^2*(r1^2+r2^2)*e3w))
672            ze3w = e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)
673            zah(ji,jj) = zah(ji,jj) + zn2 * ( wslpi(ji,jj,jk) * wslpi(ji,jj,jk)   &
674               &                            + wslpj(ji,jj,jk) * wslpj(ji,jj,jk) ) * ze3w
675            zhw(ji,jj) = zhw(ji,jj) + ze3w
676         END_3D
677      ENDIF
678
679      DO_2D( 0, 0, 0, 0 )
680         zfw = MAX( ABS( 2. * omega * SIN( rad * gphit(ji,jj) ) ) , 1.e-10 )
681         ! Rossby radius at w-point taken betwenn 2 km and  40km
682         zRo(ji,jj) = MAX(  2.e3 , MIN( .4 * zn(ji,jj) / zfw, 40.e3 )  )
683         ! Compute aeiw by multiplying Ro^2 and T^-1
684         zaeiw(ji,jj) = zRo(ji,jj) * zRo(ji,jj) * SQRT( zah(ji,jj) / zhw(ji,jj) ) * tmask(ji,jj,1)
685      END_2D
686
687      !                                         !==  Bound on eiv coeff.  ==!
688      z1_f20 = 1._wp / (  2._wp * omega * sin( rad * 20._wp )  )
689      DO_2D( 0, 0, 0, 0 )
690         zzaei = MIN( 1._wp, ABS( ff_t(ji,jj) * z1_f20 ) ) * zaeiw(ji,jj)     ! tropical decrease
691         zaeiw(ji,jj) = MIN( zzaei , paei0 )                                  ! Max value = paei0
692      END_2D
693#if defined key_mpi3
694      CALL lbc_lnk_nc_multi( 'ldftra', zaeiw(:,:), 'W', 1.0_wp )       ! lateral boundary condition
695#else
696      CALL lbc_lnk( 'ldftra', zaeiw(:,:), 'W', 1.0_wp )       ! lateral boundary condition
697#endif
698      !               
699      DO_2D( 0, 0, 0, 0 )                       !== aei at u- and v-points  ==!
700         paeiu(ji,jj,1) = 0.5_wp * ( zaeiw(ji,jj) + zaeiw(ji+1,jj  ) ) * umask(ji,jj,1)
701         paeiv(ji,jj,1) = 0.5_wp * ( zaeiw(ji,jj) + zaeiw(ji  ,jj+1) ) * vmask(ji,jj,1)
702      END_2D
703#if defined key_mpi3
704      CALL lbc_lnk_nc_multi( 'ldftra', paeiu(:,:,1), 'U', 1.0_wp , paeiv(:,:,1), 'V', 1.0_wp )      ! lateral boundary condition
705#else
706      CALL lbc_lnk_multi( 'ldftra', paeiu(:,:,1), 'U', 1.0_wp , paeiv(:,:,1), 'V', 1.0_wp )      ! lateral boundary condition
707#endif
708
709      DO jk = 2, jpkm1                          !==  deeper values equal the surface one  ==!
710         paeiu(:,:,jk) = paeiu(:,:,1) * umask(:,:,jk)
711         paeiv(:,:,jk) = paeiv(:,:,1) * vmask(:,:,jk)
712      END DO
713     
714   END SUBROUTINE ldf_eiv
715
716
717   SUBROUTINE ldf_eiv_trp( kt, kit000, pu, pv, pw, cdtype, Kmm, Krhs )
718      !!----------------------------------------------------------------------
719      !!                  ***  ROUTINE ldf_eiv_trp  ***
720      !!
721      !! ** Purpose :   add to the input ocean transport the contribution of
722      !!              the eddy induced velocity parametrization.
723      !!
724      !! ** Method  :   The eddy induced transport is computed from a flux stream-
725      !!              function which depends on the slope of iso-neutral surfaces
726      !!              (see ldf_slp). For example, in the i-k plan :
727      !!                   psi_uw = mk(aeiu) e2u mi(wslpi)   [in m3/s]
728      !!                   Utr_eiv = - dk[psi_uw]
729      !!                   Vtr_eiv = + di[psi_uw]
730      !!                ln_ldfeiv_dia = T : output the associated streamfunction,
731      !!                                    velocity and heat transport (call ldf_eiv_dia)
732      !!
733      !! ** Action  : pu, pv increased by the eiv transport
734      !!----------------------------------------------------------------------
735      INTEGER                         , INTENT(in   ) ::   kt        ! ocean time-step index
736      INTEGER                         , INTENT(in   ) ::   kit000    ! first time step index
737      INTEGER                         , INTENT(in   ) ::   Kmm, Krhs ! ocean time level indices
738      CHARACTER(len=3)                , INTENT(in   ) ::   cdtype    ! =TRA or TRC (tracer indicator)
739      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pu      ! in : 3 ocean transport components   [m3/s]
740      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pv      ! out: 3 ocean transport components   [m3/s]
741      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pw      ! increased by the eiv                [m3/s]
742      !!
743      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
744      REAL(wp) ::   zuwk, zuwk1, zuwi, zuwi1   ! local scalars
745      REAL(wp) ::   zvwk, zvwk1, zvwj, zvwj1   !   -      -
746      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpsi_uw, zpsi_vw
747      !!----------------------------------------------------------------------
748      !
749      IF( kt == kit000 )  THEN
750         IF(lwp) WRITE(numout,*)
751         IF(lwp) WRITE(numout,*) 'ldf_eiv_trp : eddy induced advection on ', cdtype,' :'
752         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   add to velocity fields the eiv component'
753      ENDIF
754
755     
756      zpsi_uw(:,:, 1 ) = 0._wp   ;   zpsi_vw(:,:, 1 ) = 0._wp
757      zpsi_uw(:,:,jpk) = 0._wp   ;   zpsi_vw(:,:,jpk) = 0._wp
758      !
759      DO_3D( 1, 0, 1, 0, 2, jpkm1 )
760         zpsi_uw(ji,jj,jk) = - r1_4 * e2u(ji,jj) * ( wslpi(ji,jj,jk  ) + wslpi(ji+1,jj,jk) )   &
761            &                                    * ( aeiu (ji,jj,jk-1) + aeiu (ji  ,jj,jk) ) * wumask(ji,jj,jk)
762         zpsi_vw(ji,jj,jk) = - r1_4 * e1v(ji,jj) * ( wslpj(ji,jj,jk  ) + wslpj(ji,jj+1,jk) )   &
763            &                                    * ( aeiv (ji,jj,jk-1) + aeiv (ji,jj  ,jk) ) * wvmask(ji,jj,jk)
764      END_3D
765      !
766      DO_3D( 1, 0, 1, 0, 1, jpkm1 )
767         pu(ji,jj,jk) = pu(ji,jj,jk) - ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji,jj,jk+1) )
768         pv(ji,jj,jk) = pv(ji,jj,jk) - ( zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj,jk+1) )
769      END_3D
770      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
771         pw(ji,jj,jk) = pw(ji,jj,jk) + (  zpsi_uw(ji,jj,jk) - zpsi_uw(ji-1,jj  ,jk)   &
772            &                             + zpsi_vw(ji,jj,jk) - zpsi_vw(ji  ,jj-1,jk) )
773      END_3D
774      !
775      !                              ! diagnose the eddy induced velocity and associated heat transport
776      IF( ln_ldfeiv_dia .AND. cdtype == 'TRA' )   CALL ldf_eiv_dia( zpsi_uw, zpsi_vw, Kmm )
777      !
778    END SUBROUTINE ldf_eiv_trp
779
780
781   SUBROUTINE ldf_eiv_dia( psi_uw, psi_vw, Kmm )
782      !!----------------------------------------------------------------------
783      !!                  ***  ROUTINE ldf_eiv_dia  ***
784      !!
785      !! ** Purpose :   diagnose the eddy induced velocity and its associated
786      !!              vertically integrated heat transport.
787      !!
788      !! ** Method :
789      !!
790      !!----------------------------------------------------------------------
791      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   psi_uw, psi_vw   ! streamfunction   [m3/s]
792      INTEGER                         , INTENT(in   ) ::   Kmm   ! ocean time level indices
793      !
794      INTEGER  ::   ji, jj, jk    ! dummy loop indices
795      REAL(wp) ::   zztmp   ! local scalar
796      REAL(wp), DIMENSION(jpi,jpj)     ::   zw2d   ! 2D workspace
797      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zw3d   ! 3D workspace
798      !!----------------------------------------------------------------------
799      !
800!!gm I don't like this routine....   Crazy  way of doing things, not optimal at all...
801!!gm     to be redesigned....   
802      !                                                  !==  eiv stream function: output  ==!
803#if defined key_mpi3
804      CALL lbc_lnk_nc_multi( 'ldftra', psi_uw, 'U', -1.0_wp , psi_vw, 'V', -1.0_wp )
805#else
806      CALL lbc_lnk_multi( 'ldftra', psi_uw, 'U', -1.0_wp , psi_vw, 'V', -1.0_wp )
807#endif
808      !
809!!gm      CALL iom_put( "psi_eiv_uw", psi_uw )                 ! output
810!!gm      CALL iom_put( "psi_eiv_vw", psi_vw )
811      !
812      !                                                  !==  eiv velocities: calculate and output  ==!
813      !
814      zw3d(:,:,jpk) = 0._wp                                    ! bottom value always 0
815      !
816      DO jk = 1, jpkm1                                         ! e2u e3u u_eiv = -dk[psi_uw]
817         zw3d(:,:,jk) = ( psi_uw(:,:,jk+1) - psi_uw(:,:,jk) ) / ( e2u(:,:) * e3u(:,:,jk,Kmm) )
818      END DO
819      CALL iom_put( "uoce_eiv", zw3d )
820      !
821      DO jk = 1, jpkm1                                         ! e1v e3v v_eiv = -dk[psi_vw]
822         zw3d(:,:,jk) = ( psi_vw(:,:,jk+1) - psi_vw(:,:,jk) ) / ( e1v(:,:) * e3v(:,:,jk,Kmm) )
823      END DO
824      CALL iom_put( "voce_eiv", zw3d )
825      !
826      DO_3D( 0, 0, 0, 0, 1, jpkm1 )                            ! e1 e2 w_eiv = dk[psix] + dk[psix]
827         zw3d(ji,jj,jk) = (  psi_vw(ji,jj,jk) - psi_vw(ji  ,jj-1,jk)    &
828            &              + psi_uw(ji,jj,jk) - psi_uw(ji-1,jj  ,jk)  ) / e1e2t(ji,jj)
829      END_3D
830#if defined key_mpi3
831      CALL lbc_lnk_nc_multi( 'ldftra', zw3d, 'T', 1.0_wp )      ! lateral boundary condition
832#else
833      CALL lbc_lnk( 'ldftra', zw3d, 'T', 1.0_wp )      ! lateral boundary condition
834#endif
835      CALL iom_put( "woce_eiv", zw3d )
836      !
837      IF( iom_use('weiv_masstr') ) THEN   ! vertical mass transport & its square value
838         zw2d(:,:) = rho0 * e1e2t(:,:)
839         DO jk = 1, jpk
840            zw3d(:,:,jk) = zw3d(:,:,jk) * zw2d(:,:)
841         END DO
842         CALL iom_put( "weiv_masstr" , zw3d ) 
843      ENDIF
844      !
845      IF( iom_use('ueiv_masstr') ) THEN
846         zw3d(:,:,:) = 0.e0
847         DO jk = 1, jpkm1
848            zw3d(:,:,jk) = rho0 * ( psi_uw(:,:,jk+1) - psi_uw(:,:,jk) ) 
849         END DO
850         CALL iom_put( "ueiv_masstr", zw3d )                  ! mass transport in i-direction
851      ENDIF
852      !
853      zztmp = 0.5_wp * rho0 * rcp 
854      IF( iom_use('ueiv_heattr') .OR. iom_use('ueiv_heattr3d') ) THEN
855        zw2d(:,:)   = 0._wp 
856        zw3d(:,:,:) = 0._wp 
857        DO_3D( 0, 0, 0, 0, 1, jpkm1 )
858           zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_uw(ji,jj,jk+1)          - psi_uw(ji  ,jj,jk)            )   &
859              &                            * ( ts    (ji,jj,jk,jp_tem,Kmm) + ts    (ji+1,jj,jk,jp_tem,Kmm) ) 
860           zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
861        END_3D
862#if defined key_mpi3
863        CALL lbc_lnk_nc_multi( 'ldftra', zw2d, 'U', -1.0_wp )
864        CALL lbc_lnk_nc_multi( 'ldftra', zw3d, 'U', -1.0_wp )
865#else
866        CALL lbc_lnk( 'ldftra', zw2d, 'U', -1.0_wp )
867        CALL lbc_lnk( 'ldftra', zw3d, 'U', -1.0_wp )
868#endif
869        CALL iom_put( "ueiv_heattr"  , zztmp * zw2d )                  ! heat transport in i-direction
870        CALL iom_put( "ueiv_heattr3d", zztmp * zw3d )                  ! heat transport in i-direction
871      ENDIF
872      !
873      IF( iom_use('veiv_masstr') ) THEN
874         zw3d(:,:,:) = 0.e0
875         DO jk = 1, jpkm1
876            zw3d(:,:,jk) = rho0 * ( psi_vw(:,:,jk+1) - psi_vw(:,:,jk) ) 
877         END DO
878         CALL iom_put( "veiv_masstr", zw3d )                  ! mass transport in i-direction
879      ENDIF
880      !
881      zw2d(:,:)   = 0._wp 
882      zw3d(:,:,:) = 0._wp 
883      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
884         zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_vw(ji,jj,jk+1)          - psi_vw(ji,jj  ,jk)            )   &
885            &                            * ( ts    (ji,jj,jk,jp_tem,Kmm) + ts    (ji,jj+1,jk,jp_tem,Kmm) ) 
886         zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
887      END_3D
888#if defined key_mpi3
889      CALL lbc_lnk_nc_multi( 'ldftra', zw2d, 'V', -1.0_wp )
890#else
891      CALL lbc_lnk( 'ldftra', zw2d, 'V', -1.0_wp )
892#endif
893      CALL iom_put( "veiv_heattr", zztmp * zw2d )                  !  heat transport in j-direction
894      CALL iom_put( "veiv_heattr", zztmp * zw3d )                  !  heat transport in j-direction
895      !
896      IF( iom_use( 'sophteiv' ) )   CALL dia_ptr_hst( jp_tem, 'eiv', 0.5 * zw3d )
897      !
898      zztmp = 0.5_wp * 0.5
899      IF( iom_use('ueiv_salttr') .OR. iom_use('ueiv_salttr3d')) THEN
900        zw2d(:,:) = 0._wp 
901        zw3d(:,:,:) = 0._wp 
902        DO_3D( 0, 0, 0, 0, 1, jpkm1 )
903           zw3d(ji,jj,jk) = zw3d(ji,jj,jk) * ( psi_uw(ji,jj,jk+1)          - psi_uw(ji  ,jj,jk)            )   &
904              &                            * ( ts    (ji,jj,jk,jp_sal,Kmm) + ts    (ji+1,jj,jk,jp_sal,Kmm) ) 
905           zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
906        END_3D
907#if defined key_mpi3
908        CALL lbc_lnk_nc_multi( 'ldftra', zw2d, 'U', -1.0_wp )
909        CALL lbc_lnk_nc_multi( 'ldftra', zw3d, 'U', -1.0_wp )
910#else
911        CALL lbc_lnk( 'ldftra', zw2d, 'U', -1.0_wp )
912        CALL lbc_lnk( 'ldftra', zw3d, 'U', -1.0_wp )
913#endif
914        CALL iom_put( "ueiv_salttr", zztmp * zw2d )                  ! salt transport in i-direction
915        CALL iom_put( "ueiv_salttr3d", zztmp * zw3d )                ! salt transport in i-direction
916      ENDIF
917      zw2d(:,:) = 0._wp 
918      zw3d(:,:,:) = 0._wp 
919      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
920         zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_vw(ji,jj,jk+1)          - psi_vw(ji,jj  ,jk)            )   &
921            &                            * ( ts    (ji,jj,jk,jp_sal,Kmm) + ts    (ji,jj+1,jk,jp_sal,Kmm) ) 
922         zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
923      END_3D
924#if defined key_mpi3
925      CALL lbc_lnk_nc_multi( 'ldftra', zw2d, 'V', -1.0_wp )
926#else
927      CALL lbc_lnk( 'ldftra', zw2d, 'V', -1.0_wp )
928#endif
929      CALL iom_put( "veiv_salttr", zztmp * zw2d )                  !  salt transport in j-direction
930      CALL iom_put( "veiv_salttr", zztmp * zw3d )                  !  salt transport in j-direction
931      !
932      IF( iom_use( 'sopsteiv' ) ) CALL dia_ptr_hst( jp_sal, 'eiv', 0.5 * zw3d )
933      !
934      !
935   END SUBROUTINE ldf_eiv_dia
936
937   !!======================================================================
938END MODULE ldftra
Note: See TracBrowser for help on using the repository browser.