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_r13383_HPC-02_Daley_Tiling/src/OCE/LDF – NEMO

source: NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/LDF/ldftra.F90 @ 13553

Last change on this file since 13553 was 13553, checked in by hadcv, 4 years ago

Merge in trunk up to [13550]

  • Property svn:keywords set to Id
File size: 52.1 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_ldfeiv .AND. .NOT.( ln_traldf_iso .OR. ln_traldf_triad ) )                &
249           &            CALL ctl_stop( 'ln_ldfeiv=T requires iso-neutral laplacian diffusion' )
250      IF( ln_isfcav .AND. ln_traldf_triad ) &
251           &            CALL ctl_stop( ' ice shelf cavity and traldf_triad not tested' )
252           !
253      IF(  nldf_tra == np_lap_i .OR. nldf_tra == np_lap_it .OR. &
254         & nldf_tra == np_blp_i .OR. nldf_tra == np_blp_it  )   l_ldfslp = .TRUE.    ! slope of neutral surfaces required
255      !
256      IF( ln_traldf_blp .AND. ( ln_traldf_iso .OR. ln_traldf_triad) ) THEN     ! iso-neutral bilaplacian need MSC
257         IF( .NOT.ln_traldf_msc )   CALL ctl_stop( 'tra_ldf_init: iso-neutral bilaplacian requires ln_traldf_msc=.true.' )
258      ENDIF
259      !
260      IF(lwp) THEN
261         WRITE(numout,*)
262         SELECT CASE( nldf_tra )
263         CASE( np_no_ldf )   ;   WRITE(numout,*) '   ==>>>   NO lateral diffusion'
264         CASE( np_lap    )   ;   WRITE(numout,*) '   ==>>>   laplacian iso-level operator'
265         CASE( np_lap_i  )   ;   WRITE(numout,*) '   ==>>>   Rotated laplacian operator (standard)'
266         CASE( np_lap_it )   ;   WRITE(numout,*) '   ==>>>   Rotated laplacian operator (triad)'
267         CASE( np_blp    )   ;   WRITE(numout,*) '   ==>>>   bilaplacian iso-level operator'
268         CASE( np_blp_i  )   ;   WRITE(numout,*) '   ==>>>   Rotated bilaplacian operator (standard)'
269         CASE( np_blp_it )   ;   WRITE(numout,*) '   ==>>>   Rotated bilaplacian operator (triad)'
270         END SELECT
271         WRITE(numout,*)
272      ENDIF
273
274      !
275      !  Space/time variation of eddy coefficients
276      ! ===========================================
277      !
278      l_ldftra_time = .FALSE.                ! no time variation except in case defined below
279      !
280      IF( ln_traldf_OFF ) THEN               !== no explicit diffusive operator  ==!
281         !
282         IF(lwp) WRITE(numout,*) '   ==>>>   No diffusive operator selected. ahtu and ahtv are not allocated'
283         RETURN
284         !
285      ELSE                                   !==  a lateral diffusion operator is used  ==!
286         !
287         !                                         ! allocate the aht arrays
288         ALLOCATE( ahtu(jpi,jpj,jpk) , ahtv(jpi,jpj,jpk) , STAT=ierr )
289         IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_tra_init: failed to allocate arrays')
290         !
291         ahtu(:,:,jpk) = 0._wp                     ! last level always 0 
292         ahtv(:,:,jpk) = 0._wp
293         !.
294         !                                         ! value of lap/blp eddy mixing coef.
295         IF(     ln_traldf_lap ) THEN   ;   zUfac = r1_2 *rn_Ud   ;   inn = 1   ;   cl_Units = ' m2/s'   !   laplacian
296         ELSEIF( ln_traldf_blp ) THEN   ;   zUfac = r1_12*rn_Ud   ;   inn = 3   ;   cl_Units = ' m4/s'   ! bilaplacian
297         ENDIF
298         aht0    = zUfac *    rn_Ld**inn              ! mixing coefficient
299         zah_max = zUfac * (ra*rad)**inn              ! maximum reachable coefficient (value at the Equator for e1=1 degree)
300         !
301         !
302         SELECT CASE(  nn_aht_ijk_t  )             !* Specification of space-time variations of ahtu, ahtv
303         !
304         CASE(   0  )      !==  constant  ==!
305            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = constant = ', aht0, cl_Units
306            ahtu(:,:,1:jpkm1) = aht0
307            ahtv(:,:,1:jpkm1) = aht0
308            !
309         CASE(  10  )      !==  fixed profile  ==!
310            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( depth )'
311            IF(lwp) WRITE(numout,*) '           surface eddy diffusivity = constant = ', aht0, cl_Units
312            ahtu(:,:,1) = aht0                        ! constant surface value
313            ahtv(:,:,1) = aht0
314            CALL ldf_c1d( 'TRA', ahtu(:,:,1), ahtv(:,:,1), ahtu, ahtv )
315            !
316         CASE ( -20 )      !== fixed horizontal shape and magnitude read in file  ==!
317            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F(i,j) read in eddy_diffusivity.nc file'
318            CALL iom_open( 'eddy_diffusivity_2D.nc', inum )
319            CALL iom_get ( inum, jpdom_global, 'ahtu_2D', ahtu(:,:,1), cd_type = 'U', psgn = 1._wp )
320            CALL iom_get ( inum, jpdom_global, 'ahtv_2D', ahtv(:,:,1), cd_type = 'V', psgn = 1._wp )
321            CALL iom_close( inum )
322            DO jk = 2, jpkm1
323               ahtu(:,:,jk) = ahtu(:,:,1)
324               ahtv(:,:,jk) = ahtv(:,:,1)
325            END DO
326            !
327         CASE(  20  )      !== fixed horizontal shape  ==!
328            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( e1, e2 ) or F( e1^3, e2^3 ) (lap or blp case)'
329            IF(lwp) WRITE(numout,*) '           using a fixed diffusive velocity = ', rn_Ud,' m/s   and   Ld = Max(e1,e2)'
330            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
331            CALL ldf_c2d( 'TRA', zUfac      , inn        , ahtu, ahtv )    ! value proportional to scale factor^inn
332            !
333         CASE(  21  )      !==  time varying 2D field  ==!
334            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( latitude, longitude, time )'
335            IF(lwp) WRITE(numout,*) '                            = F( growth rate of baroclinic instability )'
336            IF(lwp) WRITE(numout,*) '            min value = 0.2 * aht0 (with aht0= 1/2 rn_Ud*rn_Ld)'
337            IF(lwp) WRITE(numout,*) '            max value =       aei0 (with aei0=1/2 rn_Ue*Le  increased to aht0 within 20N-20S'
338            !
339            l_ldftra_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
340            !
341            IF( ln_traldf_blp )   CALL ctl_stop( 'ldf_tra_init: aht=F( growth rate of baroc. insta .)',   &
342               &                                 '              incompatible with bilaplacian operator' )
343            !
344         CASE( -30  )      !== fixed 3D shape read in file  ==!
345            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F(i,j,k) read in eddy_diffusivity.nc file'
346            CALL iom_open( 'eddy_diffusivity_3D.nc', inum )
347            CALL iom_get ( inum, jpdom_global, 'ahtu_3D', ahtu, cd_type = 'U', psgn = 1._wp )
348            CALL iom_get ( inum, jpdom_global, 'ahtv_3D', ahtv, cd_type = 'V', psgn = 1._wp )
349            CALL iom_close( inum )
350            !
351         CASE(  30  )      !==  fixed 3D shape  ==!
352            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( latitude, longitude, depth )'
353            IF(lwp) WRITE(numout,*) '           using a fixed diffusive velocity = ', rn_Ud,' m/s   and   Ld = Max(e1,e2)'
354            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
355            CALL ldf_c2d( 'TRA', zUfac      , inn        , ahtu, ahtv )    ! surface value proportional to scale factor^inn
356            CALL ldf_c1d( 'TRA', ahtu(:,:,1), ahtv(:,:,1), ahtu, ahtv )    ! reduction with depth
357            !
358         CASE(  31  )      !==  time varying 3D field  ==!
359            IF(lwp) WRITE(numout,*) '   ==>>>   eddy diffusivity = F( latitude, longitude, depth , time )'
360            IF(lwp) WRITE(numout,*) '           proportional to the velocity : 1/2 |u|e or 1/12 |u|e^3'
361            !
362            l_ldftra_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
363            !
364         CASE DEFAULT
365            CALL ctl_stop('ldf_tra_init: wrong choice for nn_aht_ijk_t, the type of space-time variation of aht')
366         END SELECT
367         !
368         IF( .NOT.l_ldftra_time ) THEN             !* No time variation
369            IF(     ln_traldf_lap ) THEN                 !   laplacian operator (mask only)
370               ahtu(:,:,1:jpkm1) =       ahtu(:,:,1:jpkm1)   * umask(:,:,1:jpkm1)
371               ahtv(:,:,1:jpkm1) =       ahtv(:,:,1:jpkm1)   * vmask(:,:,1:jpkm1)
372            ELSEIF( ln_traldf_blp ) THEN                 ! bilaplacian operator (square root + mask)
373               ahtu(:,:,1:jpkm1) = SQRT( ahtu(:,:,1:jpkm1) ) * umask(:,:,1:jpkm1)
374               ahtv(:,:,1:jpkm1) = SQRT( ahtv(:,:,1:jpkm1) ) * vmask(:,:,1:jpkm1)
375            ENDIF
376         ENDIF
377         !
378      ENDIF
379      !
380   END SUBROUTINE ldf_tra_init
381
382
383   SUBROUTINE ldf_tra( kt, Kbb, Kmm )
384      !!----------------------------------------------------------------------
385      !!                  ***  ROUTINE ldf_tra  ***
386      !!
387      !! ** Purpose :   update at kt the tracer lateral mixing coeff. (aht and aeiv)
388      !!
389      !! ** Method  : * time varying eddy diffusivity coefficients:
390      !!
391      !!    nn_aei_ijk_t = 21    aeiu, aeiv = F(i,j,  t) = F(growth rate of baroclinic instability)
392      !!                                                   with a reduction to 0 in vicinity of the Equator
393      !!    nn_aht_ijk_t = 21    ahtu, ahtv = F(i,j,  t) = F(growth rate of baroclinic instability)
394      !!
395      !!                 = 31    ahtu, ahtv = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
396      !!                                                                     or |u|e^3/12 bilaplacian operator )
397      !!
398      !!              * time varying EIV coefficients: call to ldf_eiv routine
399      !!
400      !! ** action  :   ahtu, ahtv   update at each time step   
401      !!                aeiu, aeiv      -       -     -    -   (if ln_ldfeiv=T)
402      !!----------------------------------------------------------------------
403      INTEGER, INTENT(in) ::   kt   ! time step
404      INTEGER, INTENT(in) ::   Kbb, Kmm   ! ocean time level indices
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, Kmm )
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, Kmm )
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_2D( 1, 1, 1, 1 )
433            !!gm CAUTION : here we assume lat/lon grid in 20deg N/S band (like all ORCA cfg)
434            !!     ==>>>   The Coriolis value is identical for t- & u_points, and for v- and f-points
435            zaht = ( 1._wp -  MIN( 1._wp , ABS( ff_t(ji,jj) * z1_f20 ) ) ) * zDaht
436            zahf = ( 1._wp -  MIN( 1._wp , ABS( ff_f(ji,jj) * z1_f20 ) ) ) * zDaht
437            ahtu(ji,jj,1) = (  MAX( zaht_min, ahtu(ji,jj,1) ) + zaht  )     ! min value zaht_min
438            ahtv(ji,jj,1) = (  MAX( zaht_min, ahtv(ji,jj,1) ) + zahf  )     ! increase within 20S-20N
439         END_2D
440         DO jk = 1, jpkm1                             ! deeper value = surface value + mask for all levels
441            ahtu(:,:,jk) = ahtu(:,:,1) * umask(:,:,jk)
442            ahtv(:,:,jk) = ahtv(:,:,1) * vmask(:,:,jk)
443         END DO
444         !
445      CASE(  31  )       !==  time varying 3D field  ==!   = F( local velocity )
446         IF( ln_traldf_lap     ) THEN          !   laplacian operator |u| e /12
447            DO jk = 1, jpkm1
448               ahtu(:,:,jk) = ABS( uu(:,:,jk,Kbb) ) * e1u(:,:) * r1_12   ! n.b. uu,vv are masked
449               ahtv(:,:,jk) = ABS( vv(:,:,jk,Kbb) ) * e2v(:,:) * r1_12
450            END DO
451         ELSEIF( ln_traldf_blp ) THEN      ! bilaplacian operator      sqrt( |u| e^3 /12 ) = sqrt( |u| e /12 ) * e
452            DO jk = 1, jpkm1
453               ahtu(:,:,jk) = SQRT(  ABS( uu(:,:,jk,Kbb) ) * e1u(:,:) * r1_12  ) * e1u(:,:)
454               ahtv(:,:,jk) = SQRT(  ABS( vv(:,:,jk,Kbb) ) * e2v(:,:) * r1_12  ) * e2v(:,:)
455            END DO
456         ENDIF
457         !
458      END SELECT
459      !
460      CALL iom_put( "ahtu_2d", ahtu(:,:,1) )   ! surface u-eddy diffusivity coeff.
461      CALL iom_put( "ahtv_2d", ahtv(:,:,1) )   ! surface v-eddy diffusivity coeff.
462      CALL iom_put( "ahtu_3d", ahtu(:,:,:) )   ! 3D      u-eddy diffusivity coeff.
463      CALL iom_put( "ahtv_3d", ahtv(:,:,:) )   ! 3D      v-eddy diffusivity coeff.
464      !
465      IF( ln_ldfeiv ) THEN
466        CALL iom_put( "aeiu_2d", aeiu(:,:,1) )   ! surface u-EIV coeff.
467        CALL iom_put( "aeiv_2d", aeiv(:,:,1) )   ! surface v-EIV coeff.
468        CALL iom_put( "aeiu_3d", aeiu(:,:,:) )   ! 3D      u-EIV coeff.
469        CALL iom_put( "aeiv_3d", aeiv(:,:,:) )   ! 3D      v-EIV coeff.
470      ENDIF
471      !
472   END SUBROUTINE ldf_tra
473
474
475   SUBROUTINE ldf_eiv_init
476      !!----------------------------------------------------------------------
477      !!                  ***  ROUTINE ldf_eiv_init  ***
478      !!
479      !! ** Purpose :   initialization of the eiv coeff. from namelist choices.
480      !!
481      !! ** Method  :   the eiv diffusivity coef. specification depends on:
482      !!    nn_aei_ijk_t  =  0 => = constant
483      !!                  !
484      !!                  = 10 => = F(z) : constant with a reduction of 1/4 with depth
485      !!                  !
486      !!                  =-20 => = F(i,j)   = shape read in 'eddy_diffusivity.nc' file
487      !!                  = 20    = F(i,j)   = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
488      !!                  = 21    = F(i,j,t) = F(growth rate of baroclinic instability)
489      !!                  !
490      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_diffusivity.nc' file
491      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
492      !!
493      !! ** Action  :   aeiu , aeiv   :  initialized one for all or l_ldftra_time set to true
494      !!                l_ldfeiv_time : =T if EIV coefficients vary with time
495      !!----------------------------------------------------------------------
496      INTEGER  ::   jk                     ! dummy loop indices
497      INTEGER  ::   ierr, inum, ios, inn   ! local integer
498      REAL(wp) ::   zah_max, zUfac         !   -   scalar
499      !!
500      NAMELIST/namtra_eiv/ ln_ldfeiv   , ln_ldfeiv_dia,   &   ! eddy induced velocity (eiv)
501         &                 nn_aei_ijk_t, rn_Ue, rn_Le         ! eiv  coefficient
502      !!----------------------------------------------------------------------
503      !
504      IF(lwp) THEN                      ! control print
505         WRITE(numout,*)
506         WRITE(numout,*) 'ldf_eiv_init : eddy induced velocity parametrization'
507         WRITE(numout,*) '~~~~~~~~~~~~ '
508      ENDIF
509      !
510      READ  ( numnam_ref, namtra_eiv, IOSTAT = ios, ERR = 901)
511901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtra_eiv in reference namelist' )
512      !
513      READ  ( numnam_cfg, namtra_eiv, IOSTAT = ios, ERR = 902 )
514902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtra_eiv in configuration namelist' )
515      IF(lwm)  WRITE ( numond, namtra_eiv )
516
517      IF(lwp) THEN                      ! control print
518         WRITE(numout,*) '   Namelist namtra_eiv : '
519         WRITE(numout,*) '      Eddy Induced Velocity (eiv) param.         ln_ldfeiv     = ', ln_ldfeiv
520         WRITE(numout,*) '      eiv streamfunction & velocity diag.        ln_ldfeiv_dia = ', ln_ldfeiv_dia
521         WRITE(numout,*) '      coefficients :'
522         WRITE(numout,*) '         type of time-space variation            nn_aei_ijk_t  = ', nn_aht_ijk_t
523         WRITE(numout,*) '         lateral diffusive velocity (if cst)     rn_Ue         = ', rn_Ue, ' m/s'
524         WRITE(numout,*) '         lateral diffusive length   (if cst)     rn_Le         = ', rn_Le, ' m'
525         WRITE(numout,*)
526      ENDIF
527      !
528      l_ldfeiv_time = .FALSE.       ! no time variation except in case defined below
529      !
530      !
531      IF( .NOT.ln_ldfeiv ) THEN     !== Parametrization not used  ==!
532         !
533         IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity param is NOT used'
534         ln_ldfeiv_dia = .FALSE.
535         !
536      ELSE                          !== use the parametrization  ==!
537         !
538         IF(lwp) WRITE(numout,*) '   ==>>>   use eddy induced velocity parametrization'
539         IF(lwp) WRITE(numout,*)
540         !
541         IF( ln_traldf_blp )   CALL ctl_stop( 'ldf_eiv_init: eddy induced velocity ONLY with laplacian diffusivity' )
542         !
543         !                                != allocate the aei arrays
544         ALLOCATE( aeiu(jpi,jpj,jpk), aeiv(jpi,jpj,jpk), STAT=ierr )
545         IF( ierr /= 0 )   CALL ctl_stop('STOP', 'ldf_eiv: failed to allocate arrays')
546         !
547         !                                != Specification of space-time variations of eaiu, aeiv
548         !
549         aeiu(:,:,jpk) = 0._wp               ! last level always 0 
550         aeiv(:,:,jpk) = 0._wp
551         !                                   ! value of EIV coef. (laplacian operator)
552         zUfac = r1_2 *rn_Ue                    ! velocity factor
553         inn = 1                                ! L-exponent
554         aei0    = zUfac *    rn_Le**inn        ! mixing coefficient
555         zah_max = zUfac * (ra*rad)**inn        ! maximum reachable coefficient (value at the Equator)
556
557         SELECT CASE( nn_aei_ijk_t )         !* Specification of space-time variations
558         !
559         CASE(   0  )                        !--  constant  --!
560            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = constant = ', aei0, ' m2/s'
561            aeiu(:,:,1:jpkm1) = aei0
562            aeiv(:,:,1:jpkm1) = aei0
563            !
564         CASE(  10  )                        !--  fixed profile  --!
565            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( depth )'
566            IF(lwp) WRITE(numout,*) '           surface eddy diffusivity = constant = ', aht0, ' m2/s'
567            aeiu(:,:,1) = aei0                  ! constant surface value
568            aeiv(:,:,1) = aei0
569            CALL ldf_c1d( 'TRA', aeiu(:,:,1), aeiv(:,:,1), aeiu, aeiv )
570            !
571         CASE ( -20 )                        !--  fixed horizontal shape read in file  --!
572            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F(i,j) read in eddy_diffusivity_2D.nc file'
573            CALL iom_open ( 'eddy_induced_velocity_2D.nc', inum )
574            CALL iom_get  ( inum, jpdom_global, 'aeiu', aeiu(:,:,1), cd_type = 'U', psgn = 1._wp )
575            CALL iom_get  ( inum, jpdom_global, 'aeiv', aeiv(:,:,1), cd_type = 'V', psgn = 1._wp )
576            CALL iom_close( inum )
577            DO jk = 2, jpkm1
578               aeiu(:,:,jk) = aeiu(:,:,1)
579               aeiv(:,:,jk) = aeiv(:,:,1)
580            END DO
581            !
582         CASE(  20  )                        !--  fixed horizontal shape  --!
583            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( e1, e2 )'
584            IF(lwp) WRITE(numout,*) '           using a fixed diffusive velocity = ', rn_Ue, ' m/s   and   Le = Max(e1,e2)'
585            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, ' m2/s   for e1=1°)'
586            CALL ldf_c2d( 'TRA', zUfac      , inn        , aeiu, aeiv )    ! value proportional to scale factor^inn
587            !
588         CASE(  21  )                        !--  time varying 2D field  --!
589            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( latitude, longitude, time )'
590            IF(lwp) WRITE(numout,*) '                                       = F( growth rate of baroclinic instability )'
591            IF(lwp) WRITE(numout,*) '           maximum allowed value: aei0 = ', aei0, ' m2/s'
592            !
593            l_ldfeiv_time = .TRUE.     ! will be calculated by call to ldf_tra routine in step.F90
594            !
595         CASE( -30  )                        !-- fixed 3D shape read in file  --!
596            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F(i,j,k) read in eddy_diffusivity_3D.nc file'
597            CALL iom_open ( 'eddy_induced_velocity_3D.nc', inum )
598            CALL iom_get  ( inum, jpdom_global, 'aeiu', aeiu, cd_type = 'U', psgn = 1._wp )
599            CALL iom_get  ( inum, jpdom_global, 'aeiv', aeiv, cd_type = 'V', psgn = 1._wp )
600            CALL iom_close( inum )
601            !
602         CASE(  30  )                        !--  fixed 3D shape  --!
603            IF(lwp) WRITE(numout,*) '   ==>>>   eddy induced velocity coef. = F( latitude, longitude, depth )'
604            CALL ldf_c2d( 'TRA', zUfac      , inn        , aeiu, aeiv )    ! surface value proportional to scale factor^inn
605            CALL ldf_c1d( 'TRA', aeiu(:,:,1), aeiv(:,:,1), aeiu, aeiv )    ! reduction with depth
606            !
607         CASE DEFAULT
608            CALL ctl_stop('ldf_tra_init: wrong choice for nn_aei_ijk_t, the type of space-time variation of aei')
609         END SELECT
610         !
611         IF( .NOT.l_ldfeiv_time ) THEN             !* mask if No time variation
612            DO jk = 1, jpkm1
613               aeiu(:,:,jk) = aeiu(:,:,jk) * umask(:,:,jk)
614               ahtv(:,:,jk) = ahtv(:,:,jk) * vmask(:,:,jk)
615            END DO
616         ENDIF
617         !
618      ENDIF
619      !                   
620   END SUBROUTINE ldf_eiv_init
621
622
623   SUBROUTINE ldf_eiv( kt, paei0, paeiu, paeiv, Kmm )
624      !!----------------------------------------------------------------------
625      !!                  ***  ROUTINE ldf_eiv  ***
626      !!
627      !! ** Purpose :   Compute the eddy induced velocity coefficient from the
628      !!              growth rate of baroclinic instability.
629      !!
630      !! ** Method  :   coefficient function of the growth rate of baroclinic instability
631      !!
632      !! Reference : Treguier et al. JPO 1997   ; Held and Larichev JAS 1996
633      !!----------------------------------------------------------------------
634      INTEGER                         , INTENT(in   ) ::   kt             ! ocean time-step index
635      INTEGER                         , INTENT(in   ) ::   Kmm            ! ocean time level indices
636      REAL(wp)                        , INTENT(inout) ::   paei0          ! max value            [m2/s]
637      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   paeiu, paeiv   ! eiv coefficient      [m2/s]
638      !
639      INTEGER  ::   ji, jj, jk    ! dummy loop indices
640      REAL(wp) ::   zfw, ze3w, zn2, z1_f20, zaht, zaht_min, zzaei    ! local scalars
641      REAL(wp), DIMENSION(jpi,jpj) ::   zn, zah, zhw, zRo, zaeiw   ! 2D workspace
642      !!----------------------------------------------------------------------
643      !
644      zn (:,:) = 0._wp        ! Local initialization
645      zhw(:,:) = 5._wp
646      zah(:,:) = 0._wp
647      zRo(:,:) = 0._wp
648      !                       ! Compute lateral diffusive coefficient at T-point
649      IF( ln_traldf_triad ) THEN
650         DO_3D( 0, 0, 0, 0, 1, jpk )
651            ! Take the max of N^2 and zero then take the vertical sum
652            ! of the square root of the resulting N^2 ( required to compute
653            ! internal Rossby radius Ro = .5 * sum_jpk(N) / f
654            zn2 = MAX( rn2b(ji,jj,jk), 0._wp )
655            zn(ji,jj) = zn(ji,jj) + SQRT( zn2 ) * e3w(ji,jj,jk,Kmm)
656            ! Compute elements required for the inverse time scale of baroclinic
657            ! eddies using the isopycnal slopes calculated in ldfslp.F :
658            ! T^-1 = sqrt(m_jpk(N^2*(r1^2+r2^2)*e3w))
659            ze3w = e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)
660            zah(ji,jj) = zah(ji,jj) + zn2 * wslp2(ji,jj,jk) * ze3w
661            zhw(ji,jj) = zhw(ji,jj) + ze3w
662         END_3D
663      ELSE
664         DO_3D( 0, 0, 0, 0, 1, jpk )
665            ! Take the max of N^2 and zero then take the vertical sum
666            ! of the square root of the resulting N^2 ( required to compute
667            ! internal Rossby radius Ro = .5 * sum_jpk(N) / f
668            zn2 = MAX( rn2b(ji,jj,jk), 0._wp )
669            zn(ji,jj) = zn(ji,jj) + SQRT( zn2 ) * e3w(ji,jj,jk,Kmm)
670            ! Compute elements required for the inverse time scale of baroclinic
671            ! eddies using the isopycnal slopes calculated in ldfslp.F :
672            ! T^-1 = sqrt(m_jpk(N^2*(r1^2+r2^2)*e3w))
673            ze3w = e3w(ji,jj,jk,Kmm) * wmask(ji,jj,jk)
674            zah(ji,jj) = zah(ji,jj) + zn2 * ( wslpi(ji,jj,jk) * wslpi(ji,jj,jk)   &
675               &                            + wslpj(ji,jj,jk) * wslpj(ji,jj,jk) ) * ze3w
676            zhw(ji,jj) = zhw(ji,jj) + ze3w
677         END_3D
678      ENDIF
679
680      DO_2D( 0, 0, 0, 0 )
681         zfw = MAX( ABS( 2. * omega * SIN( rad * gphit(ji,jj) ) ) , 1.e-10 )
682         ! Rossby radius at w-point taken betwenn 2 km and  40km
683         zRo(ji,jj) = MAX(  2.e3 , MIN( .4 * zn(ji,jj) / zfw, 40.e3 )  )
684         ! Compute aeiw by multiplying Ro^2 and T^-1
685         zaeiw(ji,jj) = zRo(ji,jj) * zRo(ji,jj) * SQRT( zah(ji,jj) / zhw(ji,jj) ) * tmask(ji,jj,1)
686      END_2D
687
688      !                                         !==  Bound on eiv coeff.  ==!
689      z1_f20 = 1._wp / (  2._wp * omega * sin( rad * 20._wp )  )
690      DO_2D( 0, 0, 0, 0 )
691         zzaei = MIN( 1._wp, ABS( ff_t(ji,jj) * z1_f20 ) ) * zaeiw(ji,jj)     ! tropical decrease
692         zaeiw(ji,jj) = MIN( zzaei , paei0 )                                  ! Max value = paei0
693      END_2D
694      CALL lbc_lnk( 'ldftra', zaeiw(:,:), 'W', 1.0_wp )       ! lateral boundary condition
695      !               
696      DO_2D( 0, 0, 0, 0 )                       !== aei at u- and v-points  ==!
697         paeiu(ji,jj,1) = 0.5_wp * ( zaeiw(ji,jj) + zaeiw(ji+1,jj  ) ) * umask(ji,jj,1)
698         paeiv(ji,jj,1) = 0.5_wp * ( zaeiw(ji,jj) + zaeiw(ji  ,jj+1) ) * vmask(ji,jj,1)
699      END_2D
700      CALL lbc_lnk_multi( 'ldftra', paeiu(:,:,1), 'U', 1.0_wp , paeiv(:,:,1), 'V', 1.0_wp )      ! lateral boundary condition
701
702      DO jk = 2, jpkm1                          !==  deeper values equal the surface one  ==!
703         paeiu(:,:,jk) = paeiu(:,:,1) * umask(:,:,jk)
704         paeiv(:,:,jk) = paeiv(:,:,1) * vmask(:,:,jk)
705      END DO
706     
707   END SUBROUTINE ldf_eiv
708
709
710   SUBROUTINE ldf_eiv_trp( kt, kit000, pu, pv, pw, cdtype, Kmm, Krhs )
711      !!----------------------------------------------------------------------
712      !!                  ***  ROUTINE ldf_eiv_trp  ***
713      !!
714      !! ** Purpose :   add to the input ocean transport the contribution of
715      !!              the eddy induced velocity parametrization.
716      !!
717      !! ** Method  :   The eddy induced transport is computed from a flux stream-
718      !!              function which depends on the slope of iso-neutral surfaces
719      !!              (see ldf_slp). For example, in the i-k plan :
720      !!                   psi_uw = mk(aeiu) e2u mi(wslpi)   [in m3/s]
721      !!                   Utr_eiv = - dk[psi_uw]
722      !!                   Vtr_eiv = + di[psi_uw]
723      !!                ln_ldfeiv_dia = T : output the associated streamfunction,
724      !!                                    velocity and heat transport (call ldf_eiv_dia)
725      !!
726      !! ** Action  : pu, pv increased by the eiv transport
727      !!----------------------------------------------------------------------
728      INTEGER                     , INTENT(in   ) ::   kt        ! ocean time-step index
729      INTEGER                     , INTENT(in   ) ::   kit000    ! first time step index
730      INTEGER                     , INTENT(in   ) ::   Kmm, Krhs ! ocean time level indices
731      CHARACTER(len=3)            , INTENT(in   ) ::   cdtype    ! =TRA or TRC (tracer indicator)
732      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   pu        ! in : 3 ocean transport components   [m3/s]
733      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   pv        ! out: 3 ocean transport components   [m3/s]
734      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   pw        ! increased by the eiv                [m3/s]
735      !!
736      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
737      REAL(wp) ::   zuwk, zuwk1, zuwi, zuwi1   ! local scalars
738      REAL(wp) ::   zvwk, zvwk1, zvwj, zvwj1   !   -      -
739      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk) ::   zpsi_uw, zpsi_vw
740      !!----------------------------------------------------------------------
741      !
742      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile
743         IF( kt == kit000 )  THEN
744            IF(lwp) WRITE(numout,*)
745            IF(lwp) WRITE(numout,*) 'ldf_eiv_trp : eddy induced advection on ', cdtype,' :'
746            IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   add to velocity fields the eiv component'
747         ENDIF
748      ENDIF
749
750     
751      zpsi_uw(:,:, 1 ) = 0._wp   ;   zpsi_vw(:,:, 1 ) = 0._wp
752      zpsi_uw(:,:,jpk) = 0._wp   ;   zpsi_vw(:,:,jpk) = 0._wp
753      !
754      DO_3D( 1, 0, 1, 0, 2, jpkm1 )
755         zpsi_uw(ji,jj,jk) = - r1_4 * e2u(ji,jj) * ( wslpi(ji,jj,jk  ) + wslpi(ji+1,jj,jk) )   &
756            &                                    * ( aeiu (ji,jj,jk-1) + aeiu (ji  ,jj,jk) ) * wumask(ji,jj,jk)
757         zpsi_vw(ji,jj,jk) = - r1_4 * e1v(ji,jj) * ( wslpj(ji,jj,jk  ) + wslpj(ji,jj+1,jk) )   &
758            &                                    * ( aeiv (ji,jj,jk-1) + aeiv (ji,jj  ,jk) ) * wvmask(ji,jj,jk)
759      END_3D
760      !
761      DO_3D( 1, 0, 1, 0, 1, jpkm1 )
762         pu(ji,jj,jk) = pu(ji,jj,jk) - ( zpsi_uw(ji,jj,jk) - zpsi_uw(ji,jj,jk+1) )
763         pv(ji,jj,jk) = pv(ji,jj,jk) - ( zpsi_vw(ji,jj,jk) - zpsi_vw(ji,jj,jk+1) )
764      END_3D
765      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
766         pw(ji,jj,jk) = pw(ji,jj,jk) + (  zpsi_uw(ji,jj,jk) - zpsi_uw(ji-1,jj  ,jk)   &
767            &                             + zpsi_vw(ji,jj,jk) - zpsi_vw(ji  ,jj-1,jk) )
768      END_3D
769      !
770      !                              ! diagnose the eddy induced velocity and associated heat transport
771      IF( ln_ldfeiv_dia .AND. cdtype == 'TRA' )   CALL ldf_eiv_dia( zpsi_uw, zpsi_vw, Kmm )
772      !
773    END SUBROUTINE ldf_eiv_trp
774
775
776   SUBROUTINE ldf_eiv_dia( psi_uw, psi_vw, Kmm )
777      !!----------------------------------------------------------------------
778      !!                  ***  ROUTINE ldf_eiv_dia  ***
779      !!
780      !! ** Purpose :   diagnose the eddy induced velocity and its associated
781      !!              vertically integrated heat transport.
782      !!
783      !! ** Method :
784      !!
785      !!----------------------------------------------------------------------
786      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk), INTENT(inout) ::   psi_uw, psi_vw   ! streamfunction   [m3/s]
787      INTEGER                     , INTENT(in   ) ::   Kmm   ! ocean time level indices
788      !
789      INTEGER  ::   ji, jj, jk    ! dummy loop indices
790      REAL(wp) ::   zztmp   ! local scalar
791      REAL(wp), DIMENSION(ST_2D(nn_hls))     ::   zw2d   ! 2D workspace
792      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk) ::   zw3d   ! 3D workspace
793      !!----------------------------------------------------------------------
794      !
795!!gm I don't like this routine....   Crazy  way of doing things, not optimal at all...
796!!gm     to be redesigned....   
797      !                                                  !==  eiv stream function: output  ==!
798!!gm      CALL iom_put( "psi_eiv_uw", psi_uw )                 ! output
799!!gm      CALL iom_put( "psi_eiv_vw", psi_vw )
800      !
801      !                                                  !==  eiv velocities: calculate and output  ==!
802      !
803      zw3d(:,:,jpk) = 0._wp                                    ! bottom value always 0
804      !
805      DO_3D( 0, 0, 0, 0, 1, jpkm1 )                                  ! e2u e3u u_eiv = -dk[psi_uw]
806         zw3d(ji,jj,jk) = ( psi_uw(ji,jj,jk+1) - psi_uw(ji,jj,jk) ) / ( e2u(ji,jj) * e3u(ji,jj,jk,Kmm) )
807      END_3D
808      CALL iom_put( "uoce_eiv", zw3d )
809      !
810      DO_3D( 0, 0, 0, 0, 1, jpkm1 )                                  ! e1v e3v v_eiv = -dk[psi_vw]
811         zw3d(ji,jj,jk) = ( psi_vw(ji,jj,jk+1) - psi_vw(ji,jj,jk) ) / ( e1v(ji,jj) * e3v(ji,jj,jk,Kmm) )
812      END_3D
813      CALL iom_put( "voce_eiv", zw3d )
814      !
815      DO_3D( 0, 0, 0, 0, 1, jpkm1 )                            ! e1 e2 w_eiv = dk[psix] + dk[psix]
816         zw3d(ji,jj,jk) = (  psi_vw(ji,jj,jk) - psi_vw(ji  ,jj-1,jk)    &
817            &              + psi_uw(ji,jj,jk) - psi_uw(ji-1,jj  ,jk)  ) / e1e2t(ji,jj)
818      END_3D
819      CALL iom_put( "woce_eiv", zw3d )
820      !
821      IF( iom_use('weiv_masstr') ) THEN   ! vertical mass transport & its square value
822         DO_2D( 0, 0, 0, 0 )
823            zw2d(ji,jj) = rho0 * e1e2t(ji,jj)
824         END_2D
825         DO jk = 1, jpk
826            zw3d(:,:,jk) = zw3d(:,:,jk) * zw2d(:,:)
827         END DO
828         CALL iom_put( "weiv_masstr" , zw3d ) 
829      ENDIF
830      !
831      IF( iom_use('ueiv_masstr') ) THEN
832         zw3d(:,:,:) = 0.e0
833         DO jk = 1, jpkm1
834            zw3d(:,:,jk) = rho0 * ( psi_uw(:,:,jk+1) - psi_uw(:,:,jk) ) 
835         END DO
836         CALL iom_put( "ueiv_masstr", zw3d )                  ! mass transport in i-direction
837      ENDIF
838      !
839      zztmp = 0.5_wp * rho0 * rcp 
840      IF( iom_use('ueiv_heattr') .OR. iom_use('ueiv_heattr3d') ) THEN
841        zw2d(:,:)   = 0._wp 
842        zw3d(:,:,:) = 0._wp 
843        DO_3D( 0, 0, 0, 0, 1, jpkm1 )
844           zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_uw(ji,jj,jk+1)          - psi_uw(ji  ,jj,jk)            )   &
845              &                            * ( ts    (ji,jj,jk,jp_tem,Kmm) + ts    (ji+1,jj,jk,jp_tem,Kmm) ) 
846           zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
847        END_3D
848        CALL iom_put( "ueiv_heattr"  , zztmp * zw2d )                  ! heat transport in i-direction
849        CALL iom_put( "ueiv_heattr3d", zztmp * zw3d )                  ! heat transport in i-direction
850      ENDIF
851      !
852      IF( iom_use('veiv_masstr') ) THEN
853         zw3d(:,:,:) = 0.e0
854         DO jk = 1, jpkm1
855            zw3d(:,:,jk) = rho0 * ( psi_vw(:,:,jk+1) - psi_vw(:,:,jk) ) 
856         END DO
857         CALL iom_put( "veiv_masstr", zw3d )                  ! mass transport in i-direction
858      ENDIF
859      !
860      zw2d(:,:)   = 0._wp 
861      zw3d(:,:,:) = 0._wp 
862      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
863         zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_vw(ji,jj,jk+1)          - psi_vw(ji,jj  ,jk)            )   &
864            &                            * ( ts    (ji,jj,jk,jp_tem,Kmm) + ts    (ji,jj+1,jk,jp_tem,Kmm) ) 
865         zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
866      END_3D
867      CALL iom_put( "veiv_heattr"  , zztmp * zw2d )                  !  heat transport in j-direction
868      CALL iom_put( "veiv_heattr3d", zztmp * zw3d )                  !  heat transport in j-direction
869      !
870      IF( iom_use( 'sophteiv' ) )   CALL dia_ptr_hst( jp_tem, 'eiv', 0.5 * zw3d )
871      !
872      zztmp = 0.5_wp * 0.5
873      IF( iom_use('ueiv_salttr') .OR. iom_use('ueiv_salttr3d')) THEN
874        zw2d(:,:) = 0._wp 
875        zw3d(:,:,:) = 0._wp 
876        DO_3D( 0, 0, 0, 0, 1, jpkm1 )
877           zw3d(ji,jj,jk) = zw3d(ji,jj,jk) * ( psi_uw(ji,jj,jk+1)          - psi_uw(ji  ,jj,jk)            )   &
878              &                            * ( ts    (ji,jj,jk,jp_sal,Kmm) + ts    (ji+1,jj,jk,jp_sal,Kmm) ) 
879           zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
880        END_3D
881        CALL iom_put( "ueiv_salttr", zztmp * zw2d )                  ! salt transport in i-direction
882        CALL iom_put( "ueiv_salttr3d", zztmp * zw3d )                ! salt transport in i-direction
883      ENDIF
884      zw2d(:,:) = 0._wp 
885      zw3d(:,:,:) = 0._wp 
886      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
887         zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_vw(ji,jj,jk+1)          - psi_vw(ji,jj  ,jk)            )   &
888            &                            * ( ts    (ji,jj,jk,jp_sal,Kmm) + ts    (ji,jj+1,jk,jp_sal,Kmm) ) 
889         zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk)
890      END_3D
891      CALL iom_put( "veiv_salttr"  , zztmp * zw2d )                  !  salt transport in j-direction
892      CALL iom_put( "veiv_salttr3d", zztmp * zw3d )                  !  salt transport in j-direction
893      !
894      IF( iom_use( 'sopsteiv' ) ) CALL dia_ptr_hst( jp_sal, 'eiv', 0.5 * zw3d )
895      !
896      !
897   END SUBROUTINE ldf_eiv_dia
898
899   !!======================================================================
900END MODULE ldftra
Note: See TracBrowser for help on using the repository browser.