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