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/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src/OCE/LDF – NEMO

source: NEMO/branches/2019/dev_r10973_AGRIF-01_jchanut_small_jpi_jpj/src/OCE/LDF/ldftra.F90 @ 11564

Last change on this file since 11564 was 11564, checked in by jchanut, 5 years ago

#2199, merged with trunk

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