[11305] | 1 | MODULE ablmod |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE ablmod *** |
---|
[13214] | 4 | !! Surface module : ABL computation to provide atmospheric data |
---|
[11305] | 5 | !! for surface fluxes computation |
---|
| 6 | !!====================================================================== |
---|
| 7 | !! History : 3.6 ! 2019-03 (F. Lemarié & G. Samson) Original code |
---|
| 8 | !!---------------------------------------------------------------------- |
---|
[13214] | 9 | |
---|
[11305] | 10 | !!---------------------------------------------------------------------- |
---|
| 11 | !! abl_stp : ABL single column model |
---|
| 12 | !! abl_zdf_tke : atmospheric vertical closure |
---|
| 13 | !!---------------------------------------------------------------------- |
---|
| 14 | USE abl ! ABL dynamics and tracers |
---|
| 15 | USE par_abl ! ABL constants |
---|
| 16 | |
---|
| 17 | USE phycst ! physical constants |
---|
[13214] | 18 | USE dom_oce, ONLY : tmask |
---|
[12015] | 19 | USE sbc_oce, ONLY : ght_abl, ghw_abl, e3t_abl, e3w_abl, jpka, jpkam1, rhoa |
---|
[14072] | 20 | USE sbcblk ! use rn_efac |
---|
| 21 | USE sbc_phy ! Catalog of functions for physical/meteorological parameters in the marine boundary layer |
---|
[11305] | 22 | ! |
---|
| 23 | USE prtctl ! Print control (prt_ctl routine) |
---|
| 24 | USE iom ! IOM library |
---|
| 25 | USE in_out_manager ! I/O manager |
---|
| 26 | USE lib_mpp ! MPP library |
---|
| 27 | USE timing ! Timing |
---|
| 28 | |
---|
| 29 | IMPLICIT NONE |
---|
| 30 | |
---|
| 31 | PUBLIC abl_stp ! called by sbcabl.F90 |
---|
[13214] | 32 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ustar2, zrough |
---|
[11305] | 33 | !! * Substitutions |
---|
[12340] | 34 | # include "do_loop_substitute.h90" |
---|
[11305] | 35 | |
---|
| 36 | CONTAINS |
---|
| 37 | |
---|
| 38 | |
---|
| 39 | !=================================================================================================== |
---|
[13214] | 40 | SUBROUTINE abl_stp( kt, psst, pssu, pssv, pssq, & ! in |
---|
| 41 | & pu_dta, pv_dta, pt_dta, pq_dta, & |
---|
[11334] | 42 | & pslp_dta, pgu_dta, pgv_dta, & |
---|
[14235] | 43 | & pcd_du, psen, pevp, plat, & ! in/out |
---|
[13214] | 44 | & pwndm, ptaui, ptauj, ptaum & |
---|
| 45 | #if defined key_si3 |
---|
| 46 | & , ptm_su, pssu_ice, pssv_ice & |
---|
| 47 | & , pssq_ice, pcd_du_ice, psen_ice & |
---|
| 48 | & , pevp_ice, pwndm_ice, pfrac_oce & |
---|
| 49 | & , ptaui_ice, ptauj_ice & |
---|
| 50 | #endif |
---|
| 51 | & ) |
---|
[11305] | 52 | !--------------------------------------------------------------------------------------------------- |
---|
| 53 | |
---|
| 54 | !!--------------------------------------------------------------------- |
---|
| 55 | !! *** ROUTINE abl_stp *** |
---|
| 56 | !! |
---|
[13214] | 57 | !! ** Purpose : Time-integration of the ABL model |
---|
[11305] | 58 | !! |
---|
[13214] | 59 | !! ** Method : Compute atmospheric variables : vertical turbulence |
---|
[11305] | 60 | !! + Coriolis term + newtonian relaxation |
---|
[13214] | 61 | !! |
---|
[11305] | 62 | !! ** Action : - Advance TKE to time n+1 and compute Avm_abl, Avt_abl, PBLh |
---|
| 63 | !! - Advance tracers to time n+1 (Euler backward scheme) |
---|
| 64 | !! - Compute Coriolis term with forward-backward scheme (possibly with geostrophic guide) |
---|
| 65 | !! - Advance u,v to time n+1 (Euler backward scheme) |
---|
| 66 | !! - Apply newtonian relaxation on the dynamics and the tracers |
---|
| 67 | !! - Finalize flux computation in psen, pevp, pwndm, ptaui, ptauj, ptaum |
---|
| 68 | !! |
---|
| 69 | !!---------------------------------------------------------------------- |
---|
| 70 | INTEGER , INTENT(in ) :: kt ! time step index |
---|
| 71 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: psst ! sea-surface temperature [Celsius] |
---|
| 72 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pssu ! sea-surface u (U-point) |
---|
[13214] | 73 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pssv ! sea-surface v (V-point) |
---|
[11305] | 74 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pssq ! sea-surface humidity |
---|
| 75 | REAL(wp) , INTENT(in ), DIMENSION(:,:,:) :: pu_dta ! large-scale windi |
---|
| 76 | REAL(wp) , INTENT(in ), DIMENSION(:,:,:) :: pv_dta ! large-scale windj |
---|
| 77 | REAL(wp) , INTENT(in ), DIMENSION(:,:,:) :: pgu_dta ! large-scale hpgi |
---|
| 78 | REAL(wp) , INTENT(in ), DIMENSION(:,:,:) :: pgv_dta ! large-scale hpgj |
---|
| 79 | REAL(wp) , INTENT(in ), DIMENSION(:,:,:) :: pt_dta ! large-scale pot. temp. |
---|
| 80 | REAL(wp) , INTENT(in ), DIMENSION(:,:,:) :: pq_dta ! large-scale humidity |
---|
| 81 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pslp_dta ! sea-level pressure |
---|
| 82 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pcd_du ! Cd x Du (T-point) |
---|
| 83 | REAL(wp) , INTENT(inout), DIMENSION(:,: ) :: psen ! Ch x Du |
---|
| 84 | REAL(wp) , INTENT(inout), DIMENSION(:,: ) :: pevp ! Ce x Du |
---|
[13214] | 85 | REAL(wp) , INTENT(inout), DIMENSION(:,: ) :: pwndm ! ||uwnd|| |
---|
[14235] | 86 | REAL(wp) , INTENT( out), DIMENSION(:,: ) :: plat ! latent heat flux |
---|
[11305] | 87 | REAL(wp) , INTENT( out), DIMENSION(:,: ) :: ptaui ! taux |
---|
| 88 | REAL(wp) , INTENT( out), DIMENSION(:,: ) :: ptauj ! tauy |
---|
[13214] | 89 | REAL(wp) , INTENT( out), DIMENSION(:,: ) :: ptaum ! ||tau|| |
---|
[11305] | 90 | ! |
---|
[13214] | 91 | #if defined key_si3 |
---|
[11334] | 92 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: ptm_su ! ice-surface temperature [K] |
---|
| 93 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pssu_ice ! ice-surface u (U-point) |
---|
[13214] | 94 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pssv_ice ! ice-surface v (V-point) |
---|
| 95 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pssq_ice ! ice-surface humidity |
---|
[11334] | 96 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pcd_du_ice ! Cd x Du over ice (T-point) |
---|
[11348] | 97 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: psen_ice ! Ch x Du over ice (T-point) |
---|
| 98 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pevp_ice ! Ce x Du over ice (T-point) |
---|
[13214] | 99 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pwndm_ice ! ||uwnd - uice|| |
---|
| 100 | REAL(wp) , INTENT(in ), DIMENSION(:,: ) :: pfrac_oce ! ocean fraction |
---|
[11334] | 101 | REAL(wp) , INTENT( out), DIMENSION(:,: ) :: ptaui_ice ! ice-surface taux stress (U-point) |
---|
[13214] | 102 | REAL(wp) , INTENT( out), DIMENSION(:,: ) :: ptauj_ice ! ice-surface tauy stress (V-point) |
---|
| 103 | #endif |
---|
[12015] | 104 | ! |
---|
[13214] | 105 | REAL(wp), DIMENSION(1:jpi,1:jpj ) :: zwnd_i, zwnd_j |
---|
| 106 | REAL(wp), DIMENSION(1:jpi ,2:jpka) :: zCF |
---|
[11305] | 107 | ! |
---|
[13214] | 108 | REAL(wp), DIMENSION(1:jpi ,1:jpka) :: z_elem_a |
---|
| 109 | REAL(wp), DIMENSION(1:jpi ,1:jpka) :: z_elem_b |
---|
| 110 | REAL(wp), DIMENSION(1:jpi ,1:jpka) :: z_elem_c |
---|
[11305] | 111 | ! |
---|
| 112 | INTEGER :: ji, jj, jk, jtra, jbak ! dummy loop indices |
---|
[11334] | 113 | REAL(wp) :: zztmp, zcff, ztemp, zhumi, zcff1, zztmp1, zztmp2 |
---|
| 114 | REAL(wp) :: zcff2, zfcor, zmsk, zsig, zcffu, zcffv, zzice,zzoce |
---|
[13214] | 115 | LOGICAL :: SemiImp_Cor = .TRUE. |
---|
[11305] | 116 | ! |
---|
[13214] | 117 | !!--------------------------------------------------------------------- |
---|
[11305] | 118 | ! |
---|
| 119 | IF(lwp .AND. kt == nit000) THEN ! control print |
---|
| 120 | WRITE(numout,*) |
---|
| 121 | WRITE(numout,*) 'abl_stp : ABL time stepping' |
---|
| 122 | WRITE(numout,*) '~~~~~~' |
---|
[13214] | 123 | ENDIF |
---|
[11305] | 124 | ! |
---|
| 125 | IF( kt == nit000 ) ALLOCATE ( ustar2( 1:jpi, 1:jpj ) ) |
---|
[13214] | 126 | IF( kt == nit000 ) ALLOCATE ( zrough( 1:jpi, 1:jpj ) ) |
---|
| 127 | !! Compute ustar squared as Cd || Uatm-Uoce ||^2 |
---|
| 128 | !! needed for surface boundary condition of TKE |
---|
[11305] | 129 | !! pwndm contains | U10m - U_oce | (see blk_oce_1 in sbcblk) |
---|
[13295] | 130 | DO_2D( 1, 1, 1, 1 ) |
---|
[12340] | 131 | zzoce = pCd_du (ji,jj) * pwndm (ji,jj) |
---|
[11873] | 132 | #if defined key_si3 |
---|
[13214] | 133 | zzice = pCd_du_ice(ji,jj) * pwndm_ice(ji,jj) |
---|
| 134 | ustar2(ji,jj) = zzoce * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * zzice |
---|
[11334] | 135 | #else |
---|
[13214] | 136 | ustar2(ji,jj) = zzoce |
---|
[11873] | 137 | #endif |
---|
[14072] | 138 | !#LB: sorry Cdn_oce is gone: |
---|
| 139 | !zrough(ji,jj) = ght_abl(2) * EXP( - vkarmn / SQRT( MAX( Cdn_oce(ji,jj), 1.e-4 ) ) ) !<-- recover the value of z0 from Cdn_oce |
---|
[12340] | 140 | END_2D |
---|
[14072] | 141 | |
---|
| 142 | zrough(:,:) = z0_from_Cd( ght_abl(2), pCd_du(:,:) / MAX( pwndm(:,:), 0.5_wp ) ) ! #LB: z0_from_Cd is define in sbc_phy.F90... |
---|
| 143 | |
---|
[11305] | 144 | ! |
---|
| 145 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 146 | ! ! 1 *** Advance TKE to time n+1 and compute Avm_abl, Avt_abl, PBLh |
---|
| 147 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 148 | |
---|
[13214] | 149 | CALL abl_zdf_tke( ) !--> Avm_abl, Avt_abl, pblh defined on (1,jpi) x (1,jpj) |
---|
| 150 | |
---|
[11305] | 151 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 152 | ! ! 2 *** Advance tracers to time n+1 |
---|
| 153 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[13214] | 154 | |
---|
[11305] | 155 | !------------- |
---|
| 156 | DO jj = 1, jpj ! outer loop !--> tq_abl computed on (1:jpi) x (1:jpj) |
---|
[13214] | 157 | !------------- |
---|
| 158 | ! Compute matrix elements for interior points |
---|
[11305] | 159 | DO jk = 3, jpkam1 |
---|
| 160 | DO ji = 1, jpi ! vector opt. |
---|
[13214] | 161 | z_elem_a( ji, jk ) = - rDt_abl * Avt_abl( ji, jj, jk-1 ) / e3w_abl( jk-1 ) ! lower-diagonal |
---|
| 162 | z_elem_c( ji, jk ) = - rDt_abl * Avt_abl( ji, jj, jk ) / e3w_abl( jk ) ! upper-diagonal |
---|
| 163 | z_elem_b( ji, jk ) = e3t_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk ) ! diagonal |
---|
[11305] | 164 | END DO |
---|
[13214] | 165 | END DO |
---|
| 166 | ! Boundary conditions |
---|
| 167 | DO ji = 1, jpi ! vector opt. |
---|
| 168 | ! Neumann at the bottom |
---|
| 169 | z_elem_a( ji, 2 ) = 0._wp |
---|
| 170 | z_elem_c( ji, 2 ) = - rDt_abl * Avt_abl( ji, jj, 2 ) / e3w_abl( 2 ) |
---|
[11305] | 171 | ! Homogeneous Neumann at the top |
---|
[13214] | 172 | z_elem_a( ji, jpka ) = - rDt_abl * Avt_abl( ji, jj, jpka ) / e3w_abl( jpka ) |
---|
| 173 | z_elem_c( ji, jpka ) = 0._wp |
---|
| 174 | z_elem_b( ji, jpka ) = e3t_abl( jpka ) - z_elem_a( ji, jpka ) |
---|
| 175 | END DO |
---|
[11873] | 176 | |
---|
[11305] | 177 | DO jtra = 1,jptq ! loop on active tracers |
---|
[13214] | 178 | |
---|
[11305] | 179 | DO jk = 3, jpkam1 |
---|
[13214] | 180 | !DO ji = 2, jpim1 |
---|
| 181 | DO ji = 1,jpi !!GS: to be checked if needed |
---|
| 182 | tq_abl( ji, jj, jk, nt_a, jtra ) = e3t_abl(jk) * tq_abl( ji, jj, jk, nt_n, jtra ) ! initialize right-hand-side |
---|
[11305] | 183 | END DO |
---|
| 184 | END DO |
---|
| 185 | |
---|
| 186 | IF(jtra == jp_ta) THEN |
---|
[13214] | 187 | DO ji = 1,jpi ! surface boundary condition for temperature |
---|
| 188 | zztmp1 = psen(ji, jj) |
---|
| 189 | zztmp2 = psen(ji, jj) * ( psst(ji, jj) + rt0 ) |
---|
| 190 | #if defined key_si3 |
---|
| 191 | zztmp1 = zztmp1 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * psen_ice(ji,jj) |
---|
[11334] | 192 | zztmp2 = zztmp2 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * psen_ice(ji,jj) * ptm_su(ji,jj) |
---|
[13214] | 193 | #endif |
---|
| 194 | z_elem_b( ji, 2 ) = e3t_abl( 2 ) - z_elem_c( ji, 2 ) + rDt_abl * zztmp1 |
---|
| 195 | tq_abl ( ji, jj, 2, nt_a, jtra ) = e3t_abl( 2 ) * tq_abl ( ji, jj, 2, nt_n, jtra ) + rDt_abl * zztmp2 |
---|
[11873] | 196 | tq_abl ( ji, jj, jpka, nt_a, jtra ) = e3t_abl( jpka ) * tq_abl ( ji, jj, jpka, nt_n, jtra ) |
---|
[13214] | 197 | END DO |
---|
| 198 | ELSE ! jp_qa |
---|
| 199 | DO ji = 1,jpi ! surface boundary condition for humidity |
---|
| 200 | zztmp1 = pevp(ji, jj) |
---|
| 201 | zztmp2 = pevp(ji, jj) * pssq(ji, jj) |
---|
| 202 | #if defined key_si3 |
---|
| 203 | zztmp1 = zztmp1 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pevp_ice(ji,jj) |
---|
| 204 | zztmp2 = zztmp2 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pevp_ice(ji, jj) * pssq_ice(ji, jj) |
---|
| 205 | #endif |
---|
| 206 | z_elem_b( ji, 2 ) = e3t_abl( 2 ) - z_elem_c( ji, 2 ) + rDt_abl * zztmp1 |
---|
| 207 | tq_abl ( ji, jj, 2 , nt_a, jtra ) = e3t_abl( 2 ) * tq_abl ( ji, jj, 2, nt_n, jtra ) + rDt_abl * zztmp2 |
---|
[11305] | 208 | tq_abl ( ji, jj, jpka, nt_a, jtra ) = e3t_abl( jpka ) * tq_abl ( ji, jj, jpka, nt_n, jtra ) |
---|
[13214] | 209 | END DO |
---|
[11305] | 210 | END IF |
---|
| 211 | !! |
---|
| 212 | !! Matrix inversion |
---|
| 213 | !! ---------------------------------------------------------- |
---|
[13214] | 214 | DO ji = 1,jpi |
---|
| 215 | zcff = 1._wp / z_elem_b( ji, 2 ) |
---|
| 216 | zCF ( ji, 2 ) = - zcff * z_elem_c( ji, 2 ) |
---|
| 217 | tq_abl( ji, jj, 2, nt_a, jtra ) = zcff * tq_abl( ji, jj, 2, nt_a, jtra ) |
---|
| 218 | END DO |
---|
[11305] | 219 | |
---|
[13214] | 220 | DO jk = 3, jpka |
---|
| 221 | DO ji = 1,jpi |
---|
| 222 | zcff = 1._wp / ( z_elem_b( ji, jk ) + z_elem_a( ji, jk ) * zCF( ji, jk-1 ) ) |
---|
[11305] | 223 | zCF(ji,jk) = - zcff * z_elem_c( ji, jk ) |
---|
| 224 | tq_abl(ji,jj,jk,nt_a,jtra) = zcff * ( tq_abl(ji,jj,jk ,nt_a,jtra) & |
---|
[13214] | 225 | & - z_elem_a(ji, jk) * tq_abl(ji,jj,jk-1,nt_a,jtra) ) |
---|
[11305] | 226 | END DO |
---|
| 227 | END DO |
---|
[13214] | 228 | !!FL at this point we could check positivity of tq_abl(:,:,:,nt_a,jp_qa) ... test to do ... |
---|
| 229 | DO jk = jpkam1,2,-1 |
---|
| 230 | DO ji = 1,jpi |
---|
[11305] | 231 | tq_abl(ji,jj,jk,nt_a,jtra) = tq_abl(ji,jj,jk,nt_a,jtra) + & |
---|
| 232 | & zCF(ji,jk) * tq_abl(ji,jj,jk+1,nt_a,jtra) |
---|
| 233 | END DO |
---|
| 234 | END DO |
---|
[13214] | 235 | |
---|
| 236 | END DO !<-- loop on tracers |
---|
| 237 | !! |
---|
[11305] | 238 | !------------- |
---|
[13214] | 239 | END DO ! end outer loop |
---|
| 240 | !------------- |
---|
| 241 | |
---|
[11305] | 242 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 243 | ! ! 3 *** Compute Coriolis term with geostrophic guide |
---|
[13214] | 244 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 245 | IF( SemiImp_Cor ) THEN |
---|
| 246 | |
---|
| 247 | !------------- |
---|
| 248 | DO jk = 2, jpka ! outer loop |
---|
| 249 | !------------- |
---|
| 250 | ! |
---|
| 251 | ! Advance u_abl & v_abl to time n+1 |
---|
[13295] | 252 | DO_2D( 1, 1, 1, 1 ) |
---|
[13214] | 253 | zcff = ( fft_abl(ji,jj) * rDt_abl )*( fft_abl(ji,jj) * rDt_abl ) ! (f dt)**2 |
---|
| 254 | |
---|
| 255 | u_abl( ji, jj, jk, nt_a ) = e3t_abl(jk) *( & |
---|
| 256 | & (1._wp-gamma_Cor*(1._wp-gamma_Cor)*zcff) * u_abl( ji, jj, jk, nt_n ) & |
---|
| 257 | & + rDt_abl * fft_abl(ji, jj) * v_abl( ji, jj, jk, nt_n ) ) & |
---|
| 258 | & / (1._wp + gamma_Cor*gamma_Cor*zcff) |
---|
[14072] | 259 | |
---|
[13214] | 260 | v_abl( ji, jj, jk, nt_a ) = e3t_abl(jk) *( & |
---|
| 261 | & (1._wp-gamma_Cor*(1._wp-gamma_Cor)*zcff) * v_abl( ji, jj, jk, nt_n ) & |
---|
| 262 | & - rDt_abl * fft_abl(ji, jj) * u_abl( ji, jj, jk, nt_n ) ) & |
---|
| 263 | & / (1._wp + gamma_Cor*gamma_Cor*zcff) |
---|
| 264 | END_2D |
---|
| 265 | ! |
---|
| 266 | !------------- |
---|
| 267 | END DO ! end outer loop !<-- u_abl and v_abl are properly updated on (1:jpi) x (1:jpj) |
---|
| 268 | !------------- |
---|
| 269 | ! |
---|
| 270 | IF( ln_geos_winds ) THEN |
---|
| 271 | DO jj = 1, jpj ! outer loop |
---|
| 272 | DO jk = 1, jpka |
---|
| 273 | DO ji = 1, jpi |
---|
| 274 | u_abl( ji, jj, jk, nt_a ) = u_abl( ji, jj, jk, nt_a ) & |
---|
| 275 | & - rDt_abl * e3t_abl(jk) * fft_abl(ji , jj) * pgv_dta(ji ,jj ,jk) |
---|
| 276 | v_abl( ji, jj, jk, nt_a ) = v_abl( ji, jj, jk, nt_a ) & |
---|
| 277 | & + rDt_abl * e3t_abl(jk) * fft_abl(ji, jj ) * pgu_dta(ji ,jj ,jk) |
---|
| 278 | END DO |
---|
[11305] | 279 | END DO |
---|
| 280 | END DO |
---|
[13214] | 281 | END IF |
---|
| 282 | ! |
---|
| 283 | IF( ln_hpgls_frc ) THEN |
---|
| 284 | DO jj = 1, jpj ! outer loop |
---|
| 285 | DO jk = 1, jpka |
---|
| 286 | DO ji = 1, jpi |
---|
| 287 | u_abl( ji, jj, jk, nt_a ) = u_abl( ji, jj, jk, nt_a ) - rDt_abl * e3t_abl(jk) * pgu_dta(ji,jj,jk) |
---|
| 288 | v_abl( ji, jj, jk, nt_a ) = v_abl( ji, jj, jk, nt_a ) - rDt_abl * e3t_abl(jk) * pgv_dta(ji,jj,jk) |
---|
| 289 | ENDDO |
---|
[11305] | 290 | ENDDO |
---|
| 291 | ENDDO |
---|
[13214] | 292 | END IF |
---|
| 293 | |
---|
| 294 | ELSE ! SemiImp_Cor = .FALSE. |
---|
| 295 | |
---|
| 296 | IF( ln_geos_winds ) THEN |
---|
| 297 | |
---|
| 298 | !------------- |
---|
| 299 | DO jk = 2, jpka ! outer loop |
---|
| 300 | !------------- |
---|
| 301 | ! |
---|
| 302 | IF( MOD( kt, 2 ) == 0 ) then |
---|
| 303 | ! Advance u_abl & v_abl to time n+1 |
---|
| 304 | DO jj = 1, jpj |
---|
| 305 | DO ji = 1, jpi |
---|
| 306 | zcff = fft_abl(ji,jj) * ( v_abl ( ji , jj , jk, nt_n ) - pgv_dta(ji ,jj ,jk) ) |
---|
| 307 | u_abl( ji, jj, jk, nt_a ) = u_abl( ji, jj, jk, nt_n ) + rDt_abl * zcff |
---|
| 308 | zcff = fft_abl(ji,jj) * ( u_abl ( ji , jj , jk, nt_a ) - pgu_dta(ji ,jj ,jk) ) |
---|
| 309 | v_abl( ji, jj, jk, nt_a ) = e3t_abl(jk) *( v_abl( ji, jj, jk, nt_n ) - rDt_abl * zcff ) |
---|
| 310 | u_abl( ji, jj, jk, nt_a ) = e3t_abl(jk) * u_abl( ji, jj, jk, nt_a ) |
---|
| 311 | END DO |
---|
| 312 | END DO |
---|
| 313 | ELSE |
---|
| 314 | DO jj = 1, jpj |
---|
| 315 | DO ji = 1, jpi |
---|
| 316 | zcff = fft_abl(ji,jj) * ( u_abl ( ji , jj , jk, nt_n ) - pgu_dta(ji ,jj ,jk) ) |
---|
| 317 | v_abl( ji, jj, jk, nt_a ) = v_abl( ji, jj, jk, nt_n ) - rDt_abl * zcff |
---|
| 318 | zcff = fft_abl(ji,jj) * ( v_abl ( ji , jj , jk, nt_a ) - pgv_dta(ji ,jj ,jk) ) |
---|
| 319 | u_abl( ji, jj, jk, nt_a ) = e3t_abl(jk) *( u_abl( ji, jj, jk, nt_n ) + rDt_abl * zcff ) |
---|
| 320 | v_abl( ji, jj, jk, nt_a ) = e3t_abl(jk) * v_abl( ji, jj, jk, nt_a ) |
---|
| 321 | END DO |
---|
| 322 | END DO |
---|
| 323 | END IF |
---|
| 324 | ! |
---|
| 325 | !------------- |
---|
| 326 | END DO ! end outer loop !<-- u_abl and v_abl are properly updated on (1:jpi) x (1:jpj) |
---|
| 327 | !------------- |
---|
| 328 | |
---|
| 329 | ENDIF ! ln_geos_winds |
---|
| 330 | |
---|
| 331 | ENDIF ! SemiImp_Cor |
---|
[11305] | 332 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[13214] | 333 | ! ! 4 *** Advance u,v to time n+1 |
---|
| 334 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[11305] | 335 | ! |
---|
[13214] | 336 | ! Vertical diffusion for u_abl |
---|
[11305] | 337 | !------------- |
---|
| 338 | DO jj = 1, jpj ! outer loop |
---|
[13214] | 339 | !------------- |
---|
[11305] | 340 | |
---|
| 341 | DO jk = 3, jpkam1 |
---|
[13214] | 342 | DO ji = 1, jpi |
---|
| 343 | z_elem_a( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk-1 ) / e3w_abl( jk-1 ) ! lower-diagonal |
---|
| 344 | z_elem_c( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk ) / e3w_abl( jk ) ! upper-diagonal |
---|
| 345 | z_elem_b( ji, jk ) = e3t_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk ) ! diagonal |
---|
[11305] | 346 | END DO |
---|
[13214] | 347 | END DO |
---|
| 348 | |
---|
| 349 | DO ji = 2, jpi ! boundary conditions (Avm_abl and pcd_du must be available at ji=jpi) |
---|
[11305] | 350 | !++ Surface boundary condition |
---|
[13214] | 351 | z_elem_a( ji, 2 ) = 0._wp |
---|
| 352 | z_elem_c( ji, 2 ) = - rDt_abl * Avm_abl( ji, jj, 2 ) / e3w_abl( 2 ) |
---|
[11334] | 353 | ! |
---|
[13214] | 354 | zztmp1 = pcd_du(ji, jj) |
---|
[13218] | 355 | zztmp2 = 0.5_wp * pcd_du(ji, jj) * ( pssu(ji-1, jj) + pssu(ji,jj) ) |
---|
[13214] | 356 | #if defined key_si3 |
---|
[11334] | 357 | zztmp1 = zztmp1 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) |
---|
[13218] | 358 | zzice = 0.5_wp * ( pssu_ice(ji-1, jj) + pssu_ice(ji, jj) ) |
---|
[13214] | 359 | zztmp2 = zztmp2 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) * zzice |
---|
| 360 | #endif |
---|
| 361 | z_elem_b( ji, 2 ) = e3t_abl( 2 ) - z_elem_c( ji, 2 ) + rDt_abl * zztmp1 |
---|
[12489] | 362 | u_abl( ji, jj, 2, nt_a ) = u_abl( ji, jj, 2, nt_a ) + rDt_abl * zztmp2 |
---|
[13214] | 363 | |
---|
| 364 | ! idealized test cases only |
---|
| 365 | !IF( ln_topbc_neumann ) THEN |
---|
| 366 | ! !++ Top Neumann B.C. |
---|
| 367 | ! z_elem_a( ji, jpka ) = - rDt_abl * Avm_abl( ji, jj, jpka ) / e3w_abl( jpka ) |
---|
| 368 | ! z_elem_c( ji, jpka ) = 0._wp |
---|
| 369 | ! z_elem_b( ji, jpka ) = e3t_abl( jpka ) - z_elem_a( ji, jpka ) |
---|
| 370 | ! !u_abl ( ji, jj, jpka, nt_a ) = e3t_abl( jpka ) * u_abl ( ji, jj, jpka, nt_a ) |
---|
| 371 | !ELSE |
---|
| 372 | !++ Top Dirichlet B.C. |
---|
| 373 | z_elem_a( ji, jpka ) = 0._wp |
---|
| 374 | z_elem_c( ji, jpka ) = 0._wp |
---|
| 375 | z_elem_b( ji, jpka ) = e3t_abl( jpka ) |
---|
| 376 | u_abl ( ji, jj, jpka, nt_a ) = e3t_abl( jpka ) * pu_dta(ji,jj,jk) |
---|
| 377 | !ENDIF |
---|
| 378 | |
---|
| 379 | END DO |
---|
[11305] | 380 | !! |
---|
| 381 | !! Matrix inversion |
---|
| 382 | !! ---------------------------------------------------------- |
---|
[13214] | 383 | !DO ji = 2, jpi |
---|
| 384 | DO ji = 1, jpi !!GS: TBI |
---|
[11305] | 385 | zcff = 1._wp / z_elem_b( ji, 2 ) |
---|
[13214] | 386 | zCF (ji, 2 ) = - zcff * z_elem_c( ji, 2 ) |
---|
[11305] | 387 | u_abl (ji,jj,2,nt_a) = zcff * u_abl(ji,jj,2,nt_a) |
---|
[13214] | 388 | END DO |
---|
[11305] | 389 | |
---|
[13214] | 390 | DO jk = 3, jpka |
---|
| 391 | DO ji = 2, jpi |
---|
[11305] | 392 | zcff = 1._wp / ( z_elem_b( ji, jk ) + z_elem_a( ji, jk ) * zCF (ji, jk-1 ) ) |
---|
| 393 | zCF(ji,jk) = - zcff * z_elem_c( ji, jk ) |
---|
| 394 | u_abl(ji,jj,jk,nt_a) = zcff * ( u_abl(ji,jj,jk ,nt_a) & |
---|
| 395 | & - z_elem_a(ji, jk) * u_abl(ji,jj,jk-1,nt_a) ) |
---|
| 396 | END DO |
---|
| 397 | END DO |
---|
[13214] | 398 | |
---|
| 399 | DO jk = jpkam1,2,-1 |
---|
[11305] | 400 | DO ji = 2, jpi |
---|
| 401 | u_abl(ji,jj,jk,nt_a) = u_abl(ji,jj,jk,nt_a) + zCF(ji,jk) * u_abl(ji,jj,jk+1,nt_a) |
---|
| 402 | END DO |
---|
| 403 | END DO |
---|
[13214] | 404 | |
---|
[11305] | 405 | !------------- |
---|
[13214] | 406 | END DO ! end outer loop |
---|
| 407 | !------------- |
---|
[11305] | 408 | |
---|
| 409 | ! |
---|
[13214] | 410 | ! Vertical diffusion for v_abl |
---|
[11305] | 411 | !------------- |
---|
| 412 | DO jj = 2, jpj ! outer loop |
---|
[13214] | 413 | !------------- |
---|
[11305] | 414 | ! |
---|
| 415 | DO jk = 3, jpkam1 |
---|
[13214] | 416 | DO ji = 1, jpi |
---|
| 417 | z_elem_a( ji, jk ) = -rDt_abl * Avm_abl( ji, jj, jk-1 ) / e3w_abl( jk-1 ) ! lower-diagonal |
---|
| 418 | z_elem_c( ji, jk ) = -rDt_abl * Avm_abl( ji, jj, jk ) / e3w_abl( jk ) ! upper-diagonal |
---|
| 419 | z_elem_b( ji, jk ) = e3t_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk ) ! diagonal |
---|
[11305] | 420 | END DO |
---|
| 421 | END DO |
---|
| 422 | |
---|
[13214] | 423 | DO ji = 1, jpi ! boundary conditions (Avm_abl and pcd_du must be available at jj=jpj) |
---|
[11305] | 424 | !++ Surface boundary condition |
---|
[13214] | 425 | z_elem_a( ji, 2 ) = 0._wp |
---|
| 426 | z_elem_c( ji, 2 ) = - rDt_abl * Avm_abl( ji, jj, 2 ) / e3w_abl( 2 ) |
---|
[11334] | 427 | ! |
---|
[13214] | 428 | zztmp1 = pcd_du(ji, jj) |
---|
[13218] | 429 | zztmp2 = 0.5_wp * pcd_du(ji, jj) * ( pssv(ji, jj) + pssv(ji, jj-1) ) |
---|
[13214] | 430 | #if defined key_si3 |
---|
[11334] | 431 | zztmp1 = zztmp1 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) |
---|
[13218] | 432 | zzice = 0.5_wp * ( pssv_ice(ji, jj) + pssv_ice(ji, jj-1) ) |
---|
[13214] | 433 | zztmp2 = zztmp2 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) * zzice |
---|
| 434 | #endif |
---|
| 435 | z_elem_b( ji, 2 ) = e3t_abl( 2 ) - z_elem_c( ji, 2 ) + rDt_abl * zztmp1 |
---|
[12489] | 436 | v_abl( ji, jj, 2, nt_a ) = v_abl( ji, jj, 2, nt_a ) + rDt_abl * zztmp2 |
---|
[13214] | 437 | |
---|
| 438 | ! idealized test cases only |
---|
| 439 | !IF( ln_topbc_neumann ) THEN |
---|
| 440 | ! !++ Top Neumann B.C. |
---|
| 441 | ! z_elem_a( ji, jpka ) = - rDt_abl * Avm_abl( ji, jj, jpka ) / e3w_abl( jpka ) |
---|
| 442 | ! z_elem_c( ji, jpka ) = 0._wp |
---|
| 443 | ! z_elem_b( ji, jpka ) = e3t_abl( jpka ) - z_elem_a( ji, jpka ) |
---|
| 444 | ! !v_abl ( ji, jj, jpka, nt_a ) = e3t_abl( jpka ) * v_abl ( ji, jj, jpka, nt_a ) |
---|
| 445 | !ELSE |
---|
| 446 | !++ Top Dirichlet B.C. |
---|
| 447 | z_elem_a( ji, jpka ) = 0._wp |
---|
| 448 | z_elem_c( ji, jpka ) = 0._wp |
---|
| 449 | z_elem_b( ji, jpka ) = e3t_abl( jpka ) |
---|
| 450 | v_abl ( ji, jj, jpka, nt_a ) = e3t_abl( jpka ) * pv_dta(ji,jj,jk) |
---|
| 451 | !ENDIF |
---|
| 452 | |
---|
| 453 | END DO |
---|
[11305] | 454 | !! |
---|
| 455 | !! Matrix inversion |
---|
| 456 | !! ---------------------------------------------------------- |
---|
[13214] | 457 | DO ji = 1, jpi |
---|
[11305] | 458 | zcff = 1._wp / z_elem_b( ji, 2 ) |
---|
[13214] | 459 | zCF (ji, 2 ) = - zcff * z_elem_c( ji, 2 ) |
---|
| 460 | v_abl (ji,jj,2,nt_a) = zcff * v_abl ( ji, jj, 2, nt_a ) |
---|
| 461 | END DO |
---|
[11305] | 462 | |
---|
[13214] | 463 | DO jk = 3, jpka |
---|
| 464 | DO ji = 1, jpi |
---|
[11305] | 465 | zcff = 1._wp / ( z_elem_b( ji, jk ) + z_elem_a( ji, jk ) * zCF (ji, jk-1 ) ) |
---|
| 466 | zCF(ji,jk) = - zcff * z_elem_c( ji, jk ) |
---|
| 467 | v_abl(ji,jj,jk,nt_a) = zcff * ( v_abl(ji,jj,jk ,nt_a) & |
---|
| 468 | & - z_elem_a(ji, jk) * v_abl(ji,jj,jk-1,nt_a) ) |
---|
| 469 | END DO |
---|
| 470 | END DO |
---|
[13214] | 471 | |
---|
| 472 | DO jk = jpkam1,2,-1 |
---|
| 473 | DO ji = 1, jpi |
---|
[11305] | 474 | v_abl(ji,jj,jk,nt_a) = v_abl(ji,jj,jk,nt_a) + zCF(ji,jk) * v_abl(ji,jj,jk+1,nt_a) |
---|
| 475 | END DO |
---|
| 476 | END DO |
---|
[13214] | 477 | ! |
---|
[11305] | 478 | !------------- |
---|
[13214] | 479 | END DO ! end outer loop |
---|
| 480 | !------------- |
---|
| 481 | |
---|
[11305] | 482 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[13214] | 483 | ! ! 5 *** Apply nudging on the dynamics and the tracers |
---|
| 484 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 485 | |
---|
[11305] | 486 | IF( nn_dyn_restore > 0 ) THEN |
---|
[13214] | 487 | !------------- |
---|
[11305] | 488 | DO jk = 2, jpka ! outer loop |
---|
[13214] | 489 | !------------- |
---|
[13295] | 490 | DO_2D( 0, 1, 0, 1 ) |
---|
[12340] | 491 | zcff1 = pblh( ji, jj ) |
---|
[13214] | 492 | zsig = ght_abl(jk) / MAX( jp_pblh_min, MIN( jp_pblh_max, zcff1 ) ) |
---|
| 493 | zsig = MIN( jp_bmax , MAX( zsig, jp_bmin) ) |
---|
[12340] | 494 | zmsk = msk_abl(ji,jj) |
---|
| 495 | zcff2 = jp_alp3_dyn * zsig**3 + jp_alp2_dyn * zsig**2 & |
---|
| 496 | & + jp_alp1_dyn * zsig + jp_alp0_dyn |
---|
[12489] | 497 | zcff = (1._wp-zmsk) + zmsk * zcff2 * rn_Dt ! zcff = 1 for masked points |
---|
[13214] | 498 | ! rn_Dt = rDt_abl / nn_fsbc |
---|
[12340] | 499 | zcff = zcff * rest_eq(ji,jj) |
---|
| 500 | u_abl( ji, jj, jk, nt_a ) = (1._wp - zcff ) * u_abl( ji, jj, jk, nt_a ) & |
---|
[13214] | 501 | & + zcff * pu_dta( ji, jj, jk ) |
---|
[12340] | 502 | v_abl( ji, jj, jk, nt_a ) = (1._wp - zcff ) * v_abl( ji, jj, jk, nt_a ) & |
---|
| 503 | & + zcff * pv_dta( ji, jj, jk ) |
---|
| 504 | END_2D |
---|
[11305] | 505 | !------------- |
---|
| 506 | END DO ! end outer loop |
---|
[13214] | 507 | !------------- |
---|
[11305] | 508 | END IF |
---|
| 509 | |
---|
[13214] | 510 | !------------- |
---|
[11305] | 511 | DO jk = 2, jpka ! outer loop |
---|
[13214] | 512 | !------------- |
---|
[13295] | 513 | DO_2D( 1, 1, 1, 1 ) |
---|
[12340] | 514 | zcff1 = pblh( ji, jj ) |
---|
| 515 | zsig = ght_abl(jk) / MAX( jp_pblh_min, MIN( jp_pblh_max, zcff1 ) ) |
---|
[13214] | 516 | zsig = MIN( jp_bmax , MAX( zsig, jp_bmin) ) |
---|
[12340] | 517 | zmsk = msk_abl(ji,jj) |
---|
| 518 | zcff2 = jp_alp3_tra * zsig**3 + jp_alp2_tra * zsig**2 & |
---|
| 519 | & + jp_alp1_tra * zsig + jp_alp0_tra |
---|
[12489] | 520 | zcff = (1._wp-zmsk) + zmsk * zcff2 * rn_Dt ! zcff = 1 for masked points |
---|
[13214] | 521 | ! rn_Dt = rDt_abl / nn_fsbc |
---|
[12340] | 522 | tq_abl( ji, jj, jk, nt_a, jp_ta ) = (1._wp - zcff ) * tq_abl( ji, jj, jk, nt_a, jp_ta ) & |
---|
| 523 | & + zcff * pt_dta( ji, jj, jk ) |
---|
[13214] | 524 | |
---|
[12340] | 525 | tq_abl( ji, jj, jk, nt_a, jp_qa ) = (1._wp - zcff ) * tq_abl( ji, jj, jk, nt_a, jp_qa ) & |
---|
| 526 | & + zcff * pq_dta( ji, jj, jk ) |
---|
[13214] | 527 | |
---|
[12340] | 528 | END_2D |
---|
[11305] | 529 | !------------- |
---|
| 530 | END DO ! end outer loop |
---|
[13214] | 531 | !------------- |
---|
[11305] | 532 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[14235] | 533 | ! ! 6 *** MPI exchanges & IOM outputs |
---|
[11305] | 534 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 535 | ! |
---|
[14433] | 536 | CALL lbc_lnk( 'ablmod', u_abl(:,:,:,nt_a ), 'T', -1._wp, v_abl(:,:,:,nt_a) , 'T', -1._wp ) |
---|
| 537 | CALL lbc_lnk( 'ablmod', tq_abl(:,:,:,nt_a,jp_ta), 'T', 1._wp , tq_abl(:,:,:,nt_a,jp_qa), 'T', 1._wp , kfillmode = jpfillnothing ) ! ++++ this should not be needed... |
---|
[11873] | 538 | ! |
---|
[14239] | 539 | #if defined key_xios |
---|
[13214] | 540 | ! 2D & first ABL level |
---|
| 541 | IF ( iom_use("pblh" ) ) CALL iom_put ( "pblh", pblh(:,: ) ) |
---|
[11873] | 542 | IF ( iom_use("uz1_abl") ) CALL iom_put ( "uz1_abl", u_abl(:,:,2,nt_a ) ) |
---|
| 543 | IF ( iom_use("vz1_abl") ) CALL iom_put ( "vz1_abl", v_abl(:,:,2,nt_a ) ) |
---|
| 544 | IF ( iom_use("tz1_abl") ) CALL iom_put ( "tz1_abl", tq_abl(:,:,2,nt_a,jp_ta) ) |
---|
| 545 | IF ( iom_use("qz1_abl") ) CALL iom_put ( "qz1_abl", tq_abl(:,:,2,nt_a,jp_qa) ) |
---|
| 546 | IF ( iom_use("uz1_dta") ) CALL iom_put ( "uz1_dta", pu_dta(:,:,2 ) ) |
---|
| 547 | IF ( iom_use("vz1_dta") ) CALL iom_put ( "vz1_dta", pv_dta(:,:,2 ) ) |
---|
| 548 | IF ( iom_use("tz1_dta") ) CALL iom_put ( "tz1_dta", pt_dta(:,:,2 ) ) |
---|
| 549 | IF ( iom_use("qz1_dta") ) CALL iom_put ( "qz1_dta", pq_dta(:,:,2 ) ) |
---|
[13214] | 550 | ! debug 2D |
---|
| 551 | IF( ln_geos_winds ) THEN |
---|
| 552 | IF ( iom_use("uz1_geo") ) CALL iom_put ( "uz1_geo", pgu_dta(:,:,2) ) |
---|
| 553 | IF ( iom_use("vz1_geo") ) CALL iom_put ( "vz1_geo", pgv_dta(:,:,2) ) |
---|
| 554 | END IF |
---|
| 555 | IF( ln_hpgls_frc ) THEN |
---|
| 556 | IF ( iom_use("uz1_geo") ) CALL iom_put ( "uz1_geo", pgu_dta(:,:,2)/MAX(fft_abl(:,:),2.5e-5_wp) ) |
---|
| 557 | IF ( iom_use("vz1_geo") ) CALL iom_put ( "vz1_geo", -pgv_dta(:,:,2)/MAX(fft_abl(:,:),2.5e-5_wp) ) |
---|
| 558 | END IF |
---|
| 559 | ! 3D (all ABL levels) |
---|
| 560 | IF ( iom_use("u_abl" ) ) CALL iom_put ( "u_abl" , u_abl(:,:,2:jpka,nt_a ) ) |
---|
| 561 | IF ( iom_use("v_abl" ) ) CALL iom_put ( "v_abl" , v_abl(:,:,2:jpka,nt_a ) ) |
---|
| 562 | IF ( iom_use("t_abl" ) ) CALL iom_put ( "t_abl" , tq_abl(:,:,2:jpka,nt_a,jp_ta) ) |
---|
| 563 | IF ( iom_use("q_abl" ) ) CALL iom_put ( "q_abl" , tq_abl(:,:,2:jpka,nt_a,jp_qa) ) |
---|
| 564 | IF ( iom_use("tke_abl" ) ) CALL iom_put ( "tke_abl" , tke_abl(:,:,2:jpka,nt_a ) ) |
---|
| 565 | IF ( iom_use("avm_abl" ) ) CALL iom_put ( "avm_abl" , avm_abl(:,:,2:jpka ) ) |
---|
| 566 | IF ( iom_use("avt_abl" ) ) CALL iom_put ( "avt_abl" , avt_abl(:,:,2:jpka ) ) |
---|
| 567 | IF ( iom_use("mxlm_abl") ) CALL iom_put ( "mxlm_abl", mxlm_abl(:,:,2:jpka ) ) |
---|
| 568 | IF ( iom_use("mxld_abl") ) CALL iom_put ( "mxld_abl", mxld_abl(:,:,2:jpka ) ) |
---|
| 569 | ! debug 3D |
---|
[11873] | 570 | IF ( iom_use("u_dta") ) CALL iom_put ( "u_dta", pu_dta(:,:,2:jpka) ) |
---|
| 571 | IF ( iom_use("v_dta") ) CALL iom_put ( "v_dta", pv_dta(:,:,2:jpka) ) |
---|
| 572 | IF ( iom_use("t_dta") ) CALL iom_put ( "t_dta", pt_dta(:,:,2:jpka) ) |
---|
| 573 | IF ( iom_use("q_dta") ) CALL iom_put ( "q_dta", pq_dta(:,:,2:jpka) ) |
---|
| 574 | IF( ln_geos_winds ) THEN |
---|
[13214] | 575 | IF ( iom_use("u_geo") ) CALL iom_put ( "u_geo", pgu_dta(:,:,2:jpka) ) |
---|
| 576 | IF ( iom_use("v_geo") ) CALL iom_put ( "v_geo", pgv_dta(:,:,2:jpka) ) |
---|
[11873] | 577 | END IF |
---|
| 578 | IF( ln_hpgls_frc ) THEN |
---|
[13214] | 579 | IF ( iom_use("u_geo") ) CALL iom_put ( "u_geo", pgu_dta(:,:,2:jpka)/MAX( RESHAPE( fft_abl(:,:), (/jpi,jpj,jpka-1/), fft_abl(:,:)), 2.5e-5_wp) ) |
---|
| 580 | IF ( iom_use("v_geo") ) CALL iom_put ( "v_geo", -pgv_dta(:,:,2:jpka)/MAX( RESHAPE( fft_abl(:,:), (/jpi,jpj,jpka-1/), fft_abl(:,:)), 2.5e-5_wp) ) |
---|
[11873] | 581 | END IF |
---|
[13214] | 582 | #endif |
---|
[11305] | 583 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 584 | ! ! 7 *** Finalize flux computation |
---|
[13214] | 585 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 586 | ! |
---|
[13295] | 587 | DO_2D( 1, 1, 1, 1 ) |
---|
[13214] | 588 | ztemp = tq_abl( ji, jj, 2, nt_a, jp_ta ) |
---|
| 589 | zhumi = tq_abl( ji, jj, 2, nt_a, jp_qa ) |
---|
| 590 | zcff = rho_air( ztemp, zhumi, pslp_dta( ji, jj ) ) |
---|
| 591 | psen( ji, jj ) = - cp_air(zhumi) * zcff * psen(ji,jj) * ( psst(ji,jj) + rt0 - ztemp ) !GS: negative sign to respect aerobulk convention |
---|
| 592 | pevp( ji, jj ) = rn_efac*MAX( 0._wp, zcff * pevp(ji,jj) * ( pssq(ji,jj) - zhumi ) ) |
---|
[14235] | 593 | plat( ji, jj ) = - L_vap( psst(ji,jj) ) * pevp( ji, jj ) |
---|
[13214] | 594 | rhoa( ji, jj ) = zcff |
---|
[12340] | 595 | END_2D |
---|
[13214] | 596 | |
---|
[13295] | 597 | DO_2D( 0, 1, 0, 1 ) |
---|
[14072] | 598 | zwnd_i(ji,jj) = u_abl(ji ,jj,2,nt_a) - 0.5_wp * ( pssu(ji ,jj) + pssu(ji-1,jj) ) |
---|
| 599 | zwnd_j(ji,jj) = v_abl(ji,jj ,2,nt_a) - 0.5_wp * ( pssv(ji,jj ) + pssv(ji,jj-1) ) |
---|
[12340] | 600 | END_2D |
---|
[13214] | 601 | ! |
---|
[14433] | 602 | CALL lbc_lnk( 'ablmod', zwnd_i(:,:) , 'T', -1.0_wp, zwnd_j(:,:) , 'T', -1.0_wp ) |
---|
[11305] | 603 | ! |
---|
| 604 | ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked) |
---|
[13295] | 605 | DO_2D( 1, 1, 1, 1 ) |
---|
[12340] | 606 | zcff = SQRT( zwnd_i(ji,jj) * zwnd_i(ji,jj) & |
---|
[13214] | 607 | & + zwnd_j(ji,jj) * zwnd_j(ji,jj) ) ! * msk_abl(ji,jj) |
---|
[12340] | 608 | zztmp = rhoa(ji,jj) * pcd_du(ji,jj) |
---|
[13214] | 609 | |
---|
[12340] | 610 | pwndm (ji,jj) = zcff |
---|
| 611 | ptaum (ji,jj) = zztmp * zcff |
---|
| 612 | zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj) |
---|
| 613 | zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj) |
---|
| 614 | END_2D |
---|
[11305] | 615 | ! ... utau, vtau at U- and V_points, resp. |
---|
| 616 | ! Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines |
---|
| 617 | ! Note the use of MAX(tmask(i,j),tmask(i+1,j) is to mask tau over ice shelves |
---|
[13295] | 618 | DO_2D( 0, 0, 0, 0 ) |
---|
[12340] | 619 | zcff = 0.5_wp * ( 2._wp - msk_abl(ji,jj)*msk_abl(ji+1,jj) ) |
---|
| 620 | zztmp = MAX(msk_abl(ji,jj),msk_abl(ji+1,jj)) |
---|
| 621 | ptaui(ji,jj) = zcff * zztmp * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj ) ) |
---|
| 622 | zcff = 0.5_wp * ( 2._wp - msk_abl(ji,jj)*msk_abl(ji,jj+1) ) |
---|
| 623 | zztmp = MAX(msk_abl(ji,jj),msk_abl(ji,jj+1)) |
---|
| 624 | ptauj(ji,jj) = zcff * zztmp * ( zwnd_j(ji,jj) + zwnd_j(ji ,jj+1) ) |
---|
| 625 | END_2D |
---|
[11305] | 626 | ! |
---|
[14433] | 627 | CALL lbc_lnk( 'ablmod', ptaui(:,:), 'U', -1.0_wp, ptauj(:,:), 'V', -1.0_wp ) |
---|
[11305] | 628 | |
---|
| 629 | CALL iom_put( "taum_oce", ptaum ) |
---|
| 630 | |
---|
[12236] | 631 | IF(sn_cfctl%l_prtctl) THEN |
---|
[13214] | 632 | CALL prt_ctl( tab2d_1=ptaui , clinfo1=' abl_stp: utau : ', mask1=umask, & |
---|
| 633 | & tab2d_2=ptauj , clinfo2=' vtau : ', mask2=vmask ) |
---|
| 634 | CALL prt_ctl( tab2d_1=pwndm , clinfo1=' abl_stp: wndm : ' ) |
---|
[11305] | 635 | ENDIF |
---|
[11334] | 636 | |
---|
| 637 | #if defined key_si3 |
---|
[13208] | 638 | ! ------------------------------------------------------------ ! |
---|
| 639 | ! Wind stress relative to the moving ice ( U10m - U_ice ) ! |
---|
| 640 | ! ------------------------------------------------------------ ! |
---|
[14072] | 641 | DO_2D( 0, 0, 0, 0 ) |
---|
[13208] | 642 | ptaui_ice(ji,jj) = 0.5_wp * ( rhoa(ji+1,jj) * pCd_du_ice(ji+1,jj) + rhoa(ji,jj) * pCd_du_ice(ji,jj) ) & |
---|
| 643 | & * ( 0.5_wp * ( u_abl(ji+1,jj,2,nt_a) + u_abl(ji,jj,2,nt_a) ) - pssu_ice(ji,jj) ) |
---|
| 644 | ptauj_ice(ji,jj) = 0.5_wp * ( rhoa(ji,jj+1) * pCd_du_ice(ji,jj+1) + rhoa(ji,jj) * pCd_du_ice(ji,jj) ) & |
---|
| 645 | & * ( 0.5_wp * ( v_abl(ji,jj+1,2,nt_a) + v_abl(ji,jj,2,nt_a) ) - pssv_ice(ji,jj) ) |
---|
| 646 | END_2D |
---|
[14433] | 647 | CALL lbc_lnk( 'ablmod', ptaui_ice, 'U', -1.0_wp, ptauj_ice, 'V', -1.0_wp ) |
---|
[13208] | 648 | ! |
---|
| 649 | IF(sn_cfctl%l_prtctl) CALL prt_ctl( tab2d_1=ptaui_ice , clinfo1=' abl_stp: putaui : ' & |
---|
| 650 | & , tab2d_2=ptauj_ice , clinfo2=' pvtaui : ' ) |
---|
[13214] | 651 | ! ------------------------------------------------------------ ! |
---|
| 652 | ! Wind stress relative to the moving ice ( U10m - U_ice ) ! |
---|
| 653 | ! ------------------------------------------------------------ ! |
---|
[13295] | 654 | DO_2D( 0, 0, 0, 0 ) |
---|
[13214] | 655 | |
---|
| 656 | zztmp1 = 0.5_wp * ( u_abl(ji+1,jj ,2,nt_a) + u_abl(ji,jj,2,nt_a) ) |
---|
| 657 | zztmp2 = 0.5_wp * ( v_abl(ji ,jj+1,2,nt_a) + v_abl(ji,jj,2,nt_a) ) |
---|
| 658 | |
---|
| 659 | ptaui_ice(ji,jj) = 0.5_wp * ( rhoa(ji+1,jj) * pCd_du_ice(ji+1,jj) & |
---|
| 660 | & + rhoa(ji ,jj) * pCd_du_ice(ji ,jj) ) & |
---|
[13218] | 661 | & * ( zztmp1 - pssu_ice(ji,jj) ) |
---|
[13214] | 662 | ptauj_ice(ji,jj) = 0.5_wp * ( rhoa(ji,jj+1) * pCd_du_ice(ji,jj+1) & |
---|
| 663 | & + rhoa(ji,jj ) * pCd_du_ice(ji,jj ) ) & |
---|
[13218] | 664 | & * ( zztmp2 - pssv_ice(ji,jj) ) |
---|
[13214] | 665 | END_2D |
---|
[14433] | 666 | CALL lbc_lnk( 'ablmod', ptaui_ice, 'U', -1.0_wp, ptauj_ice,'V', -1.0_wp ) |
---|
[13214] | 667 | ! |
---|
| 668 | IF(sn_cfctl%l_prtctl) THEN |
---|
| 669 | CALL prt_ctl( tab2d_1=ptaui_ice , clinfo1=' abl_stp: utau_ice : ', mask1=umask, & |
---|
| 670 | & tab2d_2=ptauj_ice , clinfo2=' vtau_ice : ', mask2=vmask ) |
---|
| 671 | END IF |
---|
[11334] | 672 | #endif |
---|
[11305] | 673 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 674 | ! ! 8 *** Swap time indices for the next timestep |
---|
[12814] | 675 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 676 | nt_n = 1 + MOD( nt_n, 2) |
---|
| 677 | nt_a = 1 + MOD( nt_a, 2) |
---|
| 678 | ! |
---|
[11305] | 679 | !--------------------------------------------------------------------------------------------------- |
---|
| 680 | END SUBROUTINE abl_stp |
---|
| 681 | !=================================================================================================== |
---|
| 682 | |
---|
| 683 | |
---|
| 684 | |
---|
| 685 | |
---|
| 686 | !=================================================================================================== |
---|
| 687 | SUBROUTINE abl_zdf_tke( ) |
---|
| 688 | !--------------------------------------------------------------------------------------------------- |
---|
| 689 | |
---|
| 690 | !!---------------------------------------------------------------------- |
---|
| 691 | !! *** ROUTINE abl_zdf_tke *** |
---|
| 692 | !! |
---|
| 693 | !! ** Purpose : Time-step Turbulente Kinetic Energy (TKE) equation |
---|
| 694 | !! |
---|
| 695 | !! ** Method : - source term due to shear |
---|
| 696 | !! - source term due to stratification |
---|
| 697 | !! - resolution of the TKE equation by inverting |
---|
| 698 | !! a tridiagonal linear system |
---|
| 699 | !! |
---|
| 700 | !! ** Action : - en : now turbulent kinetic energy) |
---|
| 701 | !! - avmu, avmv : production of TKE by shear at u and v-points |
---|
| 702 | !! (= Kz dz[Ub] * dz[Un] ) |
---|
| 703 | !! --------------------------------------------------------------------- |
---|
[13214] | 704 | INTEGER :: ji, jj, jk, tind, jbak, jkup, jkdwn |
---|
[11305] | 705 | INTEGER, DIMENSION(1:jpi ) :: ikbl |
---|
| 706 | REAL(wp) :: zcff, zcff2, ztken, zesrf, zetop, ziRic, ztv |
---|
[13214] | 707 | REAL(wp) :: zdU , zdV , zcff1, zshear, zbuoy, zsig, zustar2 |
---|
| 708 | REAL(wp) :: zdU2, zdV2, zbuoy1, zbuoy2 ! zbuoy for BL89 |
---|
| 709 | REAL(wp) :: zwndi, zwndj |
---|
[11305] | 710 | REAL(wp), DIMENSION(1:jpi, 1:jpka) :: zsh2 |
---|
| 711 | REAL(wp), DIMENSION(1:jpi,1:jpj,1:jpka) :: zbn2 |
---|
[13214] | 712 | REAL(wp), DIMENSION(1:jpi,1:jpka ) :: zFC, zRH, zCF |
---|
[11305] | 713 | REAL(wp), DIMENSION(1:jpi,1:jpka ) :: z_elem_a |
---|
| 714 | REAL(wp), DIMENSION(1:jpi,1:jpka ) :: z_elem_b |
---|
| 715 | REAL(wp), DIMENSION(1:jpi,1:jpka ) :: z_elem_c |
---|
| 716 | LOGICAL :: ln_Patankar = .FALSE. |
---|
| 717 | LOGICAL :: ln_dumpvar = .FALSE. |
---|
[13214] | 718 | LOGICAL , DIMENSION(1:jpi ) :: ln_foundl |
---|
[11305] | 719 | ! |
---|
| 720 | tind = nt_n |
---|
| 721 | ziRic = 1._wp / rn_Ric |
---|
| 722 | ! if tind = nt_a it is required to apply lbc_lnk on u_abl(nt_a) and v_abl(nt_a) |
---|
| 723 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 724 | ! ! Advance TKE equation to time n+1 |
---|
| 725 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 726 | !------------- |
---|
| 727 | DO jj = 1, jpj ! outer loop |
---|
| 728 | !------------- |
---|
| 729 | ! |
---|
[13214] | 730 | ! Compute vertical shear |
---|
[11305] | 731 | DO jk = 2, jpkam1 |
---|
[13214] | 732 | DO ji = 1, jpi |
---|
| 733 | zcff = 1.0_wp / e3w_abl( jk )**2 |
---|
| 734 | zdU = zcff* Avm_abl(ji,jj,jk) * (u_abl( ji, jj, jk+1, tind)-u_abl( ji, jj, jk , tind) )**2 |
---|
[11305] | 735 | zdV = zcff* Avm_abl(ji,jj,jk) * (v_abl( ji, jj, jk+1, tind)-v_abl( ji, jj, jk , tind) )**2 |
---|
[13214] | 736 | zsh2(ji,jk) = zdU+zdV !<-- zsh2 = Km ( ( du/dz )^2 + ( dv/dz )^2 ) |
---|
[11305] | 737 | END DO |
---|
[13214] | 738 | END DO |
---|
[11305] | 739 | ! |
---|
| 740 | ! Compute brunt-vaisala frequency |
---|
| 741 | DO jk = 2, jpkam1 |
---|
[13214] | 742 | DO ji = 1,jpi |
---|
| 743 | zcff = grav * itvref / e3w_abl( jk ) |
---|
[11305] | 744 | zcff1 = tq_abl( ji, jj, jk+1, tind, jp_ta) - tq_abl( ji, jj, jk , tind, jp_ta) |
---|
| 745 | zcff2 = tq_abl( ji, jj, jk+1, tind, jp_ta) * tq_abl( ji, jj, jk+1, tind, jp_qa) & |
---|
| 746 | & - tq_abl( ji, jj, jk , tind, jp_ta) * tq_abl( ji, jj, jk , tind, jp_qa) |
---|
| 747 | zbn2(ji,jj,jk) = zcff * ( zcff1 + rctv0 * zcff2 ) !<-- zbn2 defined on (2,jpi) |
---|
| 748 | END DO |
---|
[13214] | 749 | END DO |
---|
[11305] | 750 | ! |
---|
| 751 | ! Terms for the tridiagonal problem |
---|
| 752 | DO jk = 2, jpkam1 |
---|
[13214] | 753 | DO ji = 1, jpi |
---|
| 754 | zshear = zsh2( ji, jk ) ! zsh2 is already multiplied by Avm_abl at this point |
---|
| 755 | zsh2(ji,jk) = zsh2( ji, jk ) / Avm_abl( ji, jj, jk ) ! reformulate zsh2 as a 'true' vertical shear for PBLH computation |
---|
| 756 | zbuoy = - Avt_abl( ji, jj, jk ) * zbn2( ji, jj, jk ) |
---|
| 757 | |
---|
| 758 | z_elem_a( ji, jk ) = - 0.5_wp * rDt_abl * rn_Sch * ( Avm_abl( ji, jj, jk ) + Avm_abl( ji, jj, jk-1 ) ) / e3t_abl( jk ) ! lower-diagonal |
---|
| 759 | z_elem_c( ji, jk ) = - 0.5_wp * rDt_abl * rn_Sch * ( Avm_abl( ji, jj, jk ) + Avm_abl( ji, jj, jk+1 ) ) / e3t_abl( jk+1 ) ! upper-diagonal |
---|
[11305] | 760 | IF( (zbuoy + zshear) .gt. 0.) THEN ! Patankar trick to avoid negative values of TKE |
---|
[13214] | 761 | z_elem_b( ji, jk ) = e3w_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk ) & |
---|
| 762 | & + e3w_abl(jk) * rDt_abl * rn_Ceps * sqrt(tke_abl( ji, jj, jk, nt_n )) / mxld_abl(ji,jj,jk) ! diagonal |
---|
| 763 | tke_abl( ji, jj, jk, nt_a ) = e3w_abl(jk) * ( tke_abl( ji, jj, jk, nt_n ) + rDt_abl * ( zbuoy + zshear ) ) ! right-hand-side |
---|
[11305] | 764 | ELSE |
---|
[13214] | 765 | z_elem_b( ji, jk ) = e3w_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk ) & |
---|
| 766 | & + e3w_abl(jk) * rDt_abl * rn_Ceps * sqrt(tke_abl( ji, jj, jk, nt_n )) / mxld_abl(ji,jj,jk) & ! diagonal |
---|
| 767 | & - e3w_abl(jk) * rDt_abl * zbuoy |
---|
| 768 | tke_abl( ji, jj, jk, nt_a ) = e3w_abl(jk) * ( tke_abl( ji, jj, jk, nt_n ) + rDt_abl * zshear ) ! right-hand-side |
---|
[11305] | 769 | END IF |
---|
| 770 | END DO |
---|
[13214] | 771 | END DO |
---|
| 772 | |
---|
| 773 | DO ji = 1,jpi ! vector opt. |
---|
| 774 | zesrf = MAX( rn_Esfc * ustar2(ji,jj), tke_min ) |
---|
| 775 | zetop = tke_min |
---|
| 776 | |
---|
| 777 | z_elem_a ( ji, 1 ) = 0._wp |
---|
| 778 | z_elem_c ( ji, 1 ) = 0._wp |
---|
| 779 | z_elem_b ( ji, 1 ) = 1._wp |
---|
| 780 | tke_abl ( ji, jj, 1, nt_a ) = zesrf |
---|
| 781 | |
---|
| 782 | !++ Top Neumann B.C. |
---|
| 783 | !z_elem_a ( ji, jpka ) = - 0.5 * rDt_abl * rn_Sch * (Avm_abl(ji,jj, jpka-1 )+Avm_abl(ji,jj, jpka )) / e3t_abl( jpka ) |
---|
| 784 | !z_elem_c ( ji, jpka ) = 0._wp |
---|
| 785 | !z_elem_b ( ji, jpka ) = e3w_abl(jpka) - z_elem_a(ji, jpka ) |
---|
| 786 | !tke_abl ( ji, jj, jpka, nt_a ) = e3w_abl(jpka) * tke_abl( ji,jj, jpka, nt_n ) |
---|
| 787 | |
---|
| 788 | !++ Top Dirichlet B.C. |
---|
| 789 | z_elem_a ( ji, jpka ) = 0._wp |
---|
| 790 | z_elem_c ( ji, jpka ) = 0._wp |
---|
| 791 | z_elem_b ( ji, jpka ) = 1._wp |
---|
| 792 | tke_abl ( ji, jj, jpka, nt_a ) = zetop |
---|
| 793 | |
---|
| 794 | zbn2 ( ji, jj, 1 ) = zbn2 ( ji, jj, 2 ) |
---|
| 795 | zsh2 ( ji, 1 ) = zsh2 ( ji, 2 ) |
---|
| 796 | zbn2 ( ji, jj, jpka ) = zbn2 ( ji, jj, jpkam1 ) |
---|
| 797 | zsh2 ( ji, jpka ) = zsh2 ( ji , jpkam1 ) |
---|
| 798 | END DO |
---|
[11305] | 799 | !! |
---|
| 800 | !! Matrix inversion |
---|
| 801 | !! ---------------------------------------------------------- |
---|
| 802 | DO ji = 1,jpi |
---|
| 803 | zcff = 1._wp / z_elem_b( ji, 1 ) |
---|
[13214] | 804 | zCF (ji, 1 ) = - zcff * z_elem_c( ji, 1 ) |
---|
| 805 | tke_abl(ji,jj,1,nt_a) = zcff * tke_abl ( ji, jj, 1, nt_a ) |
---|
| 806 | END DO |
---|
[11305] | 807 | |
---|
[13214] | 808 | DO jk = 2, jpka |
---|
[11305] | 809 | DO ji = 1,jpi |
---|
| 810 | zcff = 1._wp / ( z_elem_b( ji, jk ) + z_elem_a( ji, jk ) * zCF(ji, jk-1 ) ) |
---|
| 811 | zCF(ji,jk) = - zcff * z_elem_c( ji, jk ) |
---|
| 812 | tke_abl(ji,jj,jk,nt_a) = zcff * ( tke_abl(ji,jj,jk ,nt_a) & |
---|
| 813 | & - z_elem_a(ji, jk) * tke_abl(ji,jj,jk-1,nt_a) ) |
---|
| 814 | END DO |
---|
| 815 | END DO |
---|
[13214] | 816 | |
---|
| 817 | DO jk = jpkam1,1,-1 |
---|
[11305] | 818 | DO ji = 1,jpi |
---|
| 819 | tke_abl(ji,jj,jk,nt_a) = tke_abl(ji,jj,jk,nt_a) + zCF(ji,jk) * tke_abl(ji,jj,jk+1,nt_a) |
---|
| 820 | END DO |
---|
| 821 | END DO |
---|
[13214] | 822 | |
---|
| 823 | !!FL should not be needed because of Patankar procedure |
---|
[11305] | 824 | tke_abl(2:jpi,jj,1:jpka,nt_a) = MAX( tke_abl(2:jpi,jj,1:jpka,nt_a), tke_min ) |
---|
| 825 | |
---|
| 826 | !! |
---|
| 827 | !! Diagnose PBL height |
---|
| 828 | !! ---------------------------------------------------------- |
---|
[13214] | 829 | |
---|
| 830 | |
---|
| 831 | ! |
---|
[11305] | 832 | ! arrays zRH, zFC and zCF are available at this point |
---|
| 833 | ! and zFC(:, 1 ) = 0. |
---|
| 834 | ! diagnose PBL height based on zsh2 and zbn2 |
---|
| 835 | zFC ( : ,1) = 0._wp |
---|
[13214] | 836 | ikbl( 1:jpi ) = 0 |
---|
| 837 | |
---|
[11305] | 838 | DO jk = 2,jpka |
---|
[13214] | 839 | DO ji = 1, jpi |
---|
[11305] | 840 | zcff = ghw_abl( jk-1 ) |
---|
| 841 | zcff1 = zcff / ( zcff + rn_epssfc * pblh ( ji, jj ) ) |
---|
| 842 | zcff = ghw_abl( jk ) |
---|
[11322] | 843 | zcff2 = zcff / ( zcff + rn_epssfc * pblh ( ji, jj ) ) |
---|
[11305] | 844 | zFC( ji, jk ) = zFC( ji, jk-1) + 0.5_wp * e3t_abl( jk )*( & |
---|
| 845 | zcff2 * ( zsh2( ji, jk ) - ziRic * zbn2( ji, jj, jk ) & |
---|
[11322] | 846 | - rn_Cek * ( fft_abl( ji, jj ) * fft_abl( ji, jj ) ) ) & |
---|
[11305] | 847 | + zcff1 * ( zsh2( ji, jk-1) - ziRic * zbn2( ji, jj, jk-1 ) & |
---|
[11322] | 848 | - rn_Cek * ( fft_abl( ji, jj ) * fft_abl( ji, jj ) ) ) & |
---|
[11305] | 849 | & ) |
---|
| 850 | IF( ikbl(ji) == 0 .and. zFC( ji, jk ).lt.0._wp ) ikbl(ji)=jk |
---|
| 851 | END DO |
---|
| 852 | END DO |
---|
| 853 | ! |
---|
| 854 | ! finalize the computation of the PBL height |
---|
| 855 | DO ji = 1, jpi |
---|
| 856 | jk = ikbl(ji) |
---|
| 857 | IF( jk > 2 ) THEN ! linear interpolation to get subgrid value of pblh |
---|
| 858 | pblh( ji, jj ) = ( ghw_abl( jk-1 ) * zFC( ji, jk ) & |
---|
| 859 | & - ghw_abl( jk ) * zFC( ji, jk-1 ) & |
---|
| 860 | & ) / ( zFC( ji, jk ) - zFC( ji, jk-1 ) ) |
---|
| 861 | ELSE IF( jk==2 ) THEN |
---|
| 862 | pblh( ji, jj ) = ghw_abl(2 ) |
---|
| 863 | ELSE |
---|
| 864 | pblh( ji, jj ) = ghw_abl(jpka) |
---|
[13214] | 865 | END IF |
---|
[11305] | 866 | END DO |
---|
[11322] | 867 | !------------- |
---|
[13214] | 868 | END DO |
---|
| 869 | !------------- |
---|
[11873] | 870 | ! |
---|
[13214] | 871 | ! Optional : could add pblh smoothing if pblh is noisy horizontally ... |
---|
[11873] | 872 | IF(ln_smth_pblh) THEN |
---|
[13226] | 873 | CALL lbc_lnk( 'ablmod', pblh, 'T', 1.0_wp) !, kfillmode = jpfillnothing) |
---|
[11873] | 874 | CALL smooth_pblh( pblh, msk_abl ) |
---|
[13226] | 875 | CALL lbc_lnk( 'ablmod', pblh, 'T', 1.0_wp) !, kfillmode = jpfillnothing) |
---|
[11873] | 876 | ENDIF |
---|
[11305] | 877 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 878 | ! ! Diagnostic mixing length computation |
---|
| 879 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 880 | ! |
---|
| 881 | SELECT CASE ( nn_amxl ) |
---|
| 882 | ! |
---|
[13214] | 883 | CASE ( 0 ) ! Deardroff 80 length-scale bounded by the distance to surface and bottom |
---|
| 884 | # define zlup zRH |
---|
| 885 | # define zldw zFC |
---|
[11305] | 886 | DO jj = 1, jpj ! outer loop |
---|
| 887 | ! |
---|
| 888 | DO ji = 1, jpi |
---|
[13214] | 889 | mxld_abl( ji, jj, 1 ) = mxl_min |
---|
| 890 | mxld_abl( ji, jj, jpka ) = mxl_min |
---|
| 891 | mxlm_abl( ji, jj, 1 ) = mxl_min |
---|
| 892 | mxlm_abl( ji, jj, jpka ) = mxl_min |
---|
| 893 | zldw ( ji, 1 ) = zrough(ji,jj) * rn_Lsfc |
---|
| 894 | zlup ( ji, jpka ) = mxl_min |
---|
[11305] | 895 | END DO |
---|
| 896 | ! |
---|
[13214] | 897 | DO jk = 2, jpkam1 |
---|
| 898 | DO ji = 1, jpi |
---|
| 899 | zbuoy = MAX( zbn2(ji, jj, jk), rsmall ) |
---|
| 900 | mxlm_abl( ji, jj, jk ) = MAX( mxl_min, & |
---|
| 901 | & SQRT( 2._wp * tke_abl( ji, jj, jk, nt_a ) / zbuoy ) ) |
---|
[11305] | 902 | END DO |
---|
[13214] | 903 | END DO |
---|
[11305] | 904 | ! |
---|
| 905 | ! Limit mxl |
---|
[13214] | 906 | DO jk = jpkam1,1,-1 |
---|
| 907 | DO ji = 1, jpi |
---|
| 908 | zlup(ji,jk) = MIN( zlup(ji,jk+1) + (ghw_abl(jk+1)-ghw_abl(jk)) , mxlm_abl(ji, jj, jk) ) |
---|
[11305] | 909 | END DO |
---|
[13214] | 910 | END DO |
---|
[11305] | 911 | ! |
---|
| 912 | DO jk = 2, jpka |
---|
[13214] | 913 | DO ji = 1, jpi |
---|
| 914 | zldw(ji,jk) = MIN( zldw(ji,jk-1) + (ghw_abl(jk)-ghw_abl(jk-1)) , mxlm_abl(ji, jj, jk) ) |
---|
[11305] | 915 | END DO |
---|
[13214] | 916 | END DO |
---|
[11305] | 917 | ! |
---|
[13214] | 918 | ! DO jk = 1, jpka |
---|
| 919 | ! DO ji = 1, jpi |
---|
| 920 | ! mxlm_abl( ji, jj, jk ) = SQRT( zldw( ji, jk ) * zlup( ji, jk ) ) |
---|
| 921 | ! mxld_abl( ji, jj, jk ) = MIN ( zldw( ji, jk ), zlup( ji, jk ) ) |
---|
| 922 | ! END DO |
---|
| 923 | ! END DO |
---|
| 924 | ! |
---|
[11305] | 925 | DO jk = 1, jpka |
---|
| 926 | DO ji = 1, jpi |
---|
[13214] | 927 | ! zcff = 2.*SQRT(2.)*( zldw( ji, jk )**(-2._wp/3._wp) + zlup( ji, jk )**(-2._wp/3._wp) )**(-3._wp/2._wp) |
---|
| 928 | zcff = SQRT( zldw( ji, jk ) * zlup( ji, jk ) ) |
---|
| 929 | mxlm_abl( ji, jj, jk ) = MAX( zcff, mxl_min ) |
---|
| 930 | mxld_abl( ji, jj, jk ) = MAX( MIN( zldw( ji, jk ), zlup( ji, jk ) ), mxl_min ) |
---|
[11305] | 931 | END DO |
---|
[13214] | 932 | END DO |
---|
[11305] | 933 | ! |
---|
| 934 | END DO |
---|
[13214] | 935 | # undef zlup |
---|
| 936 | # undef zldw |
---|
[11305] | 937 | ! |
---|
| 938 | ! |
---|
[13214] | 939 | CASE ( 1 ) ! Modified Deardroff 80 length-scale bounded by the distance to surface and bottom |
---|
| 940 | # define zlup zRH |
---|
| 941 | # define zldw zFC |
---|
| 942 | DO jj = 1, jpj ! outer loop |
---|
[11305] | 943 | ! |
---|
[13214] | 944 | DO jk = 2, jpkam1 |
---|
| 945 | DO ji = 1,jpi |
---|
| 946 | zcff = 1.0_wp / e3w_abl( jk )**2 |
---|
| 947 | zdU = zcff* (u_abl( ji, jj, jk+1, tind)-u_abl( ji, jj, jk , tind) )**2 |
---|
| 948 | zdV = zcff* (v_abl( ji, jj, jk+1, tind)-v_abl( ji, jj, jk , tind) )**2 |
---|
| 949 | zsh2(ji,jk) = SQRT(zdU+zdV) !<-- zsh2 = SQRT ( ( du/dz )^2 + ( dv/dz )^2 ) |
---|
| 950 | ENDDO |
---|
| 951 | ENDDO |
---|
| 952 | ! |
---|
| 953 | DO ji = 1, jpi |
---|
| 954 | zcff = zrough(ji,jj) * rn_Lsfc |
---|
| 955 | mxld_abl ( ji, jj, 1 ) = zcff |
---|
| 956 | mxld_abl ( ji, jj, jpka ) = mxl_min |
---|
| 957 | mxlm_abl ( ji, jj, 1 ) = zcff |
---|
| 958 | mxlm_abl ( ji, jj, jpka ) = mxl_min |
---|
| 959 | zldw ( ji, 1 ) = zcff |
---|
| 960 | zlup ( ji, jpka ) = mxl_min |
---|
| 961 | END DO |
---|
| 962 | ! |
---|
| 963 | DO jk = 2, jpkam1 |
---|
| 964 | DO ji = 1, jpi |
---|
| 965 | zbuoy = MAX( zbn2(ji, jj, jk), rsmall ) |
---|
[13226] | 966 | zcff = 2.0_wp*SQRT(tke_abl( ji, jj, jk, nt_a )) / ( rn_Rod*zsh2(ji,jk) & |
---|
| 967 | & + SQRT(rn_Rod*rn_Rod*zsh2(ji,jk)*zsh2(ji,jk)+2.0_wp*zbuoy ) ) |
---|
[13214] | 968 | mxlm_abl( ji, jj, jk ) = MAX( mxl_min, zcff ) |
---|
[11305] | 969 | END DO |
---|
[13214] | 970 | END DO |
---|
| 971 | ! |
---|
| 972 | ! Limit mxl |
---|
| 973 | DO jk = jpkam1,1,-1 |
---|
| 974 | DO ji = 1, jpi |
---|
| 975 | zlup(ji,jk) = MIN( zlup(ji,jk+1) + (ghw_abl(jk+1)-ghw_abl(jk)) , mxlm_abl(ji, jj, jk) ) |
---|
| 976 | END DO |
---|
| 977 | END DO |
---|
| 978 | ! |
---|
| 979 | DO jk = 2, jpka |
---|
| 980 | DO ji = 1, jpi |
---|
| 981 | zldw(ji,jk) = MIN( zldw(ji,jk-1) + (ghw_abl(jk)-ghw_abl(jk-1)) , mxlm_abl(ji, jj, jk) ) |
---|
| 982 | END DO |
---|
| 983 | END DO |
---|
| 984 | ! |
---|
| 985 | DO jk = 1, jpka |
---|
| 986 | DO ji = 1, jpi |
---|
| 987 | !mxlm_abl( ji, jj, jk ) = SQRT( zldw( ji, jk ) * zlup( ji, jk ) ) |
---|
| 988 | !zcff = 2.*SQRT(2.)*( zldw( ji, jk )**(-2._wp/3._wp) + zlup( ji, jk )**(-2._wp/3._wp) )**(-3._wp/2._wp) |
---|
| 989 | zcff = SQRT( zldw( ji, jk ) * zlup( ji, jk ) ) |
---|
| 990 | mxlm_abl( ji, jj, jk ) = MAX( zcff, mxl_min ) |
---|
| 991 | !mxld_abl( ji, jj, jk ) = MIN( zldw( ji, jk ), zlup( ji, jk ) ) |
---|
| 992 | mxld_abl( ji, jj, jk ) = MAX( MIN( zldw( ji, jk ), zlup( ji, jk ) ), mxl_min ) |
---|
| 993 | END DO |
---|
| 994 | END DO |
---|
| 995 | ! |
---|
| 996 | END DO |
---|
| 997 | # undef zlup |
---|
| 998 | # undef zldw |
---|
[11305] | 999 | ! |
---|
| 1000 | CASE ( 2 ) ! Bougeault & Lacarrere 89 length-scale |
---|
| 1001 | ! |
---|
[13214] | 1002 | # define zlup zRH |
---|
| 1003 | # define zldw zFC |
---|
[11305] | 1004 | ! zCF is used for matrix inversion |
---|
[13214] | 1005 | ! |
---|
[11305] | 1006 | DO jj = 1, jpj ! outer loop |
---|
[13214] | 1007 | |
---|
| 1008 | DO ji = 1, jpi |
---|
| 1009 | zcff = zrough(ji,jj) * rn_Lsfc |
---|
| 1010 | zlup( ji, 1 ) = zcff |
---|
| 1011 | zldw( ji, 1 ) = zcff |
---|
[11305] | 1012 | zlup( ji, jpka ) = mxl_min |
---|
[13214] | 1013 | zldw( ji, jpka ) = mxl_min |
---|
| 1014 | END DO |
---|
| 1015 | |
---|
[11305] | 1016 | DO jk = 2,jpka-1 |
---|
| 1017 | DO ji = 1, jpi |
---|
| 1018 | zlup(ji,jk) = ghw_abl(jpka) - ghw_abl(jk) |
---|
[13214] | 1019 | zldw(ji,jk) = ghw_abl(jk ) - ghw_abl( 1) |
---|
[11305] | 1020 | END DO |
---|
[13214] | 1021 | END DO |
---|
[11305] | 1022 | !! |
---|
| 1023 | !! BL89 search for lup |
---|
[13214] | 1024 | !! ---------------------------------------------------------- |
---|
| 1025 | DO jk=2,jpka-1 |
---|
[11305] | 1026 | ! |
---|
| 1027 | DO ji = 1, jpi |
---|
| 1028 | zCF(ji,1:jpka) = 0._wp |
---|
| 1029 | zCF(ji, jk ) = - tke_abl( ji, jj, jk, nt_a ) |
---|
| 1030 | ln_foundl(ji ) = .false. |
---|
[13214] | 1031 | END DO |
---|
| 1032 | ! |
---|
[11305] | 1033 | DO jkup=jk+1,jpka-1 |
---|
| 1034 | DO ji = 1, jpi |
---|
[13214] | 1035 | zbuoy1 = MAX( zbn2(ji,jj,jkup ), rsmall ) |
---|
| 1036 | zbuoy2 = MAX( zbn2(ji,jj,jkup-1), rsmall ) |
---|
[11305] | 1037 | zCF (ji,jkup) = zCF (ji,jkup-1) + 0.5_wp * e3t_abl(jkup) * & |
---|
[13214] | 1038 | & ( zbuoy1*(ghw_abl(jkup )-ghw_abl(jk)) & |
---|
| 1039 | & + zbuoy2*(ghw_abl(jkup-1)-ghw_abl(jk)) ) |
---|
[11305] | 1040 | IF( zCF (ji,jkup) * zCF (ji,jkup-1) .le. 0._wp .and. .not. ln_foundl(ji) ) THEN |
---|
| 1041 | zcff2 = ghw_abl(jkup ) - ghw_abl(jk) |
---|
[13214] | 1042 | zcff1 = ghw_abl(jkup-1) - ghw_abl(jk) |
---|
[11305] | 1043 | zcff = ( zcff1 * zCF(ji,jkup) - zcff2 * zCF(ji,jkup-1) ) / & |
---|
[13214] | 1044 | & ( zCF(ji,jkup) - zCF(ji,jkup-1) ) |
---|
| 1045 | zlup(ji,jk) = zcff |
---|
| 1046 | zlup(ji,jk) = ghw_abl(jkup ) - ghw_abl(jk) |
---|
[11305] | 1047 | ln_foundl(ji) = .true. |
---|
| 1048 | END IF |
---|
| 1049 | END DO |
---|
| 1050 | END DO |
---|
| 1051 | ! |
---|
[13214] | 1052 | END DO |
---|
[11305] | 1053 | !! |
---|
| 1054 | !! BL89 search for ldwn |
---|
[13214] | 1055 | !! ---------------------------------------------------------- |
---|
| 1056 | DO jk=2,jpka-1 |
---|
[11305] | 1057 | ! |
---|
| 1058 | DO ji = 1, jpi |
---|
| 1059 | zCF(ji,1:jpka) = 0._wp |
---|
| 1060 | zCF(ji, jk ) = - tke_abl( ji, jj, jk, nt_a ) |
---|
| 1061 | ln_foundl(ji ) = .false. |
---|
[13214] | 1062 | END DO |
---|
| 1063 | ! |
---|
[11305] | 1064 | DO jkdwn=jk-1,1,-1 |
---|
[13214] | 1065 | DO ji = 1, jpi |
---|
| 1066 | zbuoy1 = MAX( zbn2(ji,jj,jkdwn+1), rsmall ) |
---|
| 1067 | zbuoy2 = MAX( zbn2(ji,jj,jkdwn ), rsmall ) |
---|
[11305] | 1068 | zCF (ji,jkdwn) = zCF (ji,jkdwn+1) + 0.5_wp * e3t_abl(jkdwn+1) & |
---|
[13214] | 1069 | & * ( zbuoy1*(ghw_abl(jk)-ghw_abl(jkdwn+1)) & |
---|
| 1070 | + zbuoy2*(ghw_abl(jk)-ghw_abl(jkdwn )) ) |
---|
| 1071 | IF(zCF (ji,jkdwn) * zCF (ji,jkdwn+1) .le. 0._wp .and. .not. ln_foundl(ji) ) THEN |
---|
[11305] | 1072 | zcff2 = ghw_abl(jk) - ghw_abl(jkdwn+1) |
---|
[13214] | 1073 | zcff1 = ghw_abl(jk) - ghw_abl(jkdwn ) |
---|
[11305] | 1074 | zcff = ( zcff1 * zCF(ji,jkdwn+1) - zcff2 * zCF(ji,jkdwn) ) / & |
---|
[13214] | 1075 | & ( zCF(ji,jkdwn+1) - zCF(ji,jkdwn) ) |
---|
| 1076 | zldw(ji,jk) = zcff |
---|
| 1077 | zldw(ji,jk) = ghw_abl(jk) - ghw_abl(jkdwn ) |
---|
| 1078 | ln_foundl(ji) = .true. |
---|
| 1079 | END IF |
---|
[11305] | 1080 | END DO |
---|
| 1081 | END DO |
---|
[13214] | 1082 | ! |
---|
[11305] | 1083 | END DO |
---|
| 1084 | |
---|
| 1085 | DO jk = 1, jpka |
---|
[13214] | 1086 | DO ji = 1, jpi |
---|
| 1087 | !zcff = 2.*SQRT(2.)*( zldw( ji, jk )**(-2._wp/3._wp) + zlup( ji, jk )**(-2._wp/3._wp) )**(-3._wp/2._wp) |
---|
| 1088 | zcff = SQRT( zldw( ji, jk ) * zlup( ji, jk ) ) |
---|
| 1089 | mxlm_abl( ji, jj, jk ) = MAX( zcff, mxl_min ) |
---|
| 1090 | mxld_abl( ji, jj, jk ) = MAX( MIN( zldw( ji, jk ), zlup( ji, jk ) ), mxl_min ) |
---|
[11305] | 1091 | END DO |
---|
[13214] | 1092 | END DO |
---|
[11305] | 1093 | |
---|
| 1094 | END DO |
---|
[13214] | 1095 | # undef zlup |
---|
| 1096 | # undef zldw |
---|
[11305] | 1097 | ! |
---|
[13214] | 1098 | CASE ( 3 ) ! Bougeault & Lacarrere 89 length-scale |
---|
| 1099 | ! |
---|
| 1100 | # define zlup zRH |
---|
| 1101 | # define zldw zFC |
---|
| 1102 | ! zCF is used for matrix inversion |
---|
| 1103 | ! |
---|
| 1104 | DO jj = 1, jpj ! outer loop |
---|
| 1105 | ! |
---|
| 1106 | DO jk = 2, jpkam1 |
---|
| 1107 | DO ji = 1,jpi |
---|
| 1108 | zcff = 1.0_wp / e3w_abl( jk )**2 |
---|
| 1109 | zdU = zcff* (u_abl( ji, jj, jk+1, tind)-u_abl( ji, jj, jk , tind) )**2 |
---|
| 1110 | zdV = zcff* (v_abl( ji, jj, jk+1, tind)-v_abl( ji, jj, jk , tind) )**2 |
---|
| 1111 | zsh2(ji,jk) = SQRT(zdU+zdV) !<-- zsh2 = SQRT ( ( du/dz )^2 + ( dv/dz )^2 ) |
---|
| 1112 | ENDDO |
---|
| 1113 | ENDDO |
---|
| 1114 | zsh2(:, 1) = zsh2( :, 2) |
---|
| 1115 | zsh2(:, jpka) = zsh2( :, jpkam1) |
---|
| 1116 | |
---|
| 1117 | DO ji = 1, jpi |
---|
| 1118 | zcff = zrough(ji,jj) * rn_Lsfc |
---|
| 1119 | zlup( ji, 1 ) = zcff |
---|
| 1120 | zldw( ji, 1 ) = zcff |
---|
| 1121 | zlup( ji, jpka ) = mxl_min |
---|
| 1122 | zldw( ji, jpka ) = mxl_min |
---|
| 1123 | END DO |
---|
| 1124 | |
---|
| 1125 | DO jk = 2,jpka-1 |
---|
| 1126 | DO ji = 1, jpi |
---|
| 1127 | zlup(ji,jk) = ghw_abl(jpka) - ghw_abl(jk) |
---|
| 1128 | zldw(ji,jk) = ghw_abl(jk ) - ghw_abl( 1) |
---|
| 1129 | END DO |
---|
| 1130 | END DO |
---|
| 1131 | !! |
---|
| 1132 | !! BL89 search for lup |
---|
| 1133 | !! ---------------------------------------------------------- |
---|
| 1134 | DO jk=2,jpka-1 |
---|
| 1135 | ! |
---|
| 1136 | DO ji = 1, jpi |
---|
| 1137 | zCF(ji,1:jpka) = 0._wp |
---|
| 1138 | zCF(ji, jk ) = - tke_abl( ji, jj, jk, nt_a ) |
---|
| 1139 | ln_foundl(ji ) = .false. |
---|
| 1140 | END DO |
---|
| 1141 | ! |
---|
| 1142 | DO jkup=jk+1,jpka-1 |
---|
| 1143 | DO ji = 1, jpi |
---|
| 1144 | zbuoy1 = MAX( zbn2(ji,jj,jkup ), rsmall ) |
---|
| 1145 | zbuoy2 = MAX( zbn2(ji,jj,jkup-1), rsmall ) |
---|
| 1146 | zCF (ji,jkup) = zCF (ji,jkup-1) + 0.5_wp * e3t_abl(jkup) * & |
---|
| 1147 | & ( zbuoy1*(ghw_abl(jkup )-ghw_abl(jk)) & |
---|
| 1148 | & + zbuoy2*(ghw_abl(jkup-1)-ghw_abl(jk)) ) & |
---|
| 1149 | & + 0.5_wp * e3t_abl(jkup) * rn_Rod * & |
---|
| 1150 | & ( SQRT(tke_abl( ji, jj, jkup , nt_a ))*zsh2(ji,jkup ) & |
---|
| 1151 | & + SQRT(tke_abl( ji, jj, jkup-1, nt_a ))*zsh2(ji,jkup-1) ) |
---|
| 1152 | |
---|
| 1153 | IF( zCF (ji,jkup) * zCF (ji,jkup-1) .le. 0._wp .and. .not. ln_foundl(ji) ) THEN |
---|
| 1154 | zcff2 = ghw_abl(jkup ) - ghw_abl(jk) |
---|
| 1155 | zcff1 = ghw_abl(jkup-1) - ghw_abl(jk) |
---|
| 1156 | zcff = ( zcff1 * zCF(ji,jkup) - zcff2 * zCF(ji,jkup-1) ) / & |
---|
| 1157 | & ( zCF(ji,jkup) - zCF(ji,jkup-1) ) |
---|
| 1158 | zlup(ji,jk) = zcff |
---|
| 1159 | zlup(ji,jk) = ghw_abl(jkup ) - ghw_abl(jk) |
---|
| 1160 | ln_foundl(ji) = .true. |
---|
| 1161 | END IF |
---|
| 1162 | END DO |
---|
| 1163 | END DO |
---|
| 1164 | ! |
---|
| 1165 | END DO |
---|
| 1166 | !! |
---|
| 1167 | !! BL89 search for ldwn |
---|
| 1168 | !! ---------------------------------------------------------- |
---|
| 1169 | DO jk=2,jpka-1 |
---|
| 1170 | ! |
---|
| 1171 | DO ji = 1, jpi |
---|
| 1172 | zCF(ji,1:jpka) = 0._wp |
---|
| 1173 | zCF(ji, jk ) = - tke_abl( ji, jj, jk, nt_a ) |
---|
| 1174 | ln_foundl(ji ) = .false. |
---|
| 1175 | END DO |
---|
| 1176 | ! |
---|
| 1177 | DO jkdwn=jk-1,1,-1 |
---|
| 1178 | DO ji = 1, jpi |
---|
| 1179 | zbuoy1 = MAX( zbn2(ji,jj,jkdwn+1), rsmall ) |
---|
| 1180 | zbuoy2 = MAX( zbn2(ji,jj,jkdwn ), rsmall ) |
---|
| 1181 | zCF (ji,jkdwn) = zCF (ji,jkdwn+1) + 0.5_wp * e3t_abl(jkdwn+1) & |
---|
| 1182 | & * (zbuoy1*(ghw_abl(jk)-ghw_abl(jkdwn+1)) & |
---|
| 1183 | & +zbuoy2*(ghw_abl(jk)-ghw_abl(jkdwn )) ) & |
---|
| 1184 | & + 0.5_wp * e3t_abl(jkup) * rn_Rod * & |
---|
| 1185 | & ( SQRT(tke_abl( ji, jj, jkdwn+1, nt_a ))*zsh2(ji,jkdwn+1) & |
---|
| 1186 | & + SQRT(tke_abl( ji, jj, jkdwn , nt_a ))*zsh2(ji,jkdwn ) ) |
---|
| 1187 | |
---|
| 1188 | IF(zCF (ji,jkdwn) * zCF (ji,jkdwn+1) .le. 0._wp .and. .not. ln_foundl(ji) ) THEN |
---|
| 1189 | zcff2 = ghw_abl(jk) - ghw_abl(jkdwn+1) |
---|
| 1190 | zcff1 = ghw_abl(jk) - ghw_abl(jkdwn ) |
---|
| 1191 | zcff = ( zcff1 * zCF(ji,jkdwn+1) - zcff2 * zCF(ji,jkdwn) ) / & |
---|
| 1192 | & ( zCF(ji,jkdwn+1) - zCF(ji,jkdwn) ) |
---|
| 1193 | zldw(ji,jk) = zcff |
---|
| 1194 | zldw(ji,jk) = ghw_abl(jk) - ghw_abl(jkdwn ) |
---|
| 1195 | ln_foundl(ji) = .true. |
---|
| 1196 | END IF |
---|
| 1197 | END DO |
---|
| 1198 | END DO |
---|
| 1199 | ! |
---|
| 1200 | END DO |
---|
| 1201 | |
---|
| 1202 | DO jk = 1, jpka |
---|
| 1203 | DO ji = 1, jpi |
---|
| 1204 | !zcff = 2.*SQRT(2.)*( zldw( ji, jk )**(-2._wp/3._wp) + zlup( ji, jk )**(-2._wp/3._wp) )**(-3._wp/2._wp) |
---|
| 1205 | zcff = SQRT( zldw( ji, jk ) * zlup( ji, jk ) ) |
---|
| 1206 | mxlm_abl( ji, jj, jk ) = MAX( zcff, mxl_min ) |
---|
| 1207 | mxld_abl( ji, jj, jk ) = MAX( MIN( zldw( ji, jk ), zlup( ji, jk ) ), mxl_min ) |
---|
| 1208 | END DO |
---|
| 1209 | END DO |
---|
| 1210 | |
---|
| 1211 | END DO |
---|
| 1212 | # undef zlup |
---|
| 1213 | # undef zldw |
---|
| 1214 | ! |
---|
| 1215 | ! |
---|
| 1216 | END SELECT |
---|
[11305] | 1217 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 1218 | ! ! Finalize the computation of turbulent visc./diff. |
---|
| 1219 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[13214] | 1220 | |
---|
[11305] | 1221 | !------------- |
---|
| 1222 | DO jj = 1, jpj ! outer loop |
---|
| 1223 | !------------- |
---|
[13214] | 1224 | DO jk = 1, jpka |
---|
[11305] | 1225 | DO ji = 1, jpi ! vector opt. |
---|
[13214] | 1226 | zcff = MAX( rn_phimax, rn_Ric * mxlm_abl( ji, jj, jk ) * mxld_abl( ji, jj, jk ) & |
---|
| 1227 | & * MAX( zbn2(ji, jj, jk), rsmall ) / tke_abl( ji, jj, jk, nt_a ) ) |
---|
| 1228 | zcff2 = 1. / ( 1. + zcff ) !<-- phi_z(z) |
---|
| 1229 | zcff = mxlm_abl( ji, jj, jk ) * SQRT( tke_abl( ji, jj, jk, nt_a ) ) |
---|
| 1230 | !!FL: MAX function probably useless because of the definition of mxl_min |
---|
[11305] | 1231 | Avm_abl( ji, jj, jk ) = MAX( rn_Cm * zcff , avm_bak ) |
---|
[13214] | 1232 | Avt_abl( ji, jj, jk ) = MAX( rn_Ct * zcff * zcff2 , avt_bak ) |
---|
[11305] | 1233 | END DO |
---|
| 1234 | END DO |
---|
| 1235 | !------------- |
---|
[13214] | 1236 | END DO |
---|
| 1237 | !------------- |
---|
[11305] | 1238 | |
---|
| 1239 | !--------------------------------------------------------------------------------------------------- |
---|
| 1240 | END SUBROUTINE abl_zdf_tke |
---|
| 1241 | !=================================================================================================== |
---|
| 1242 | |
---|
| 1243 | |
---|
[11322] | 1244 | !=================================================================================================== |
---|
| 1245 | SUBROUTINE smooth_pblh( pvar2d, msk ) |
---|
| 1246 | !--------------------------------------------------------------------------------------------------- |
---|
| 1247 | |
---|
| 1248 | !!---------------------------------------------------------------------- |
---|
| 1249 | !! *** ROUTINE smooth_pblh *** |
---|
| 1250 | !! |
---|
| 1251 | !! ** Purpose : 2D Hanning filter on atmospheric PBL height |
---|
| 1252 | !! |
---|
| 1253 | !! --------------------------------------------------------------------- |
---|
[13214] | 1254 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in ) :: msk |
---|
| 1255 | REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: pvar2d |
---|
[11322] | 1256 | INTEGER :: ji,jj |
---|
[13214] | 1257 | REAL(wp) :: smth_a, smth_b |
---|
| 1258 | REAL(wp), DIMENSION(jpi,jpj) :: zdX,zdY,zFX,zFY |
---|
| 1259 | REAL(wp) :: zumsk,zvmsk |
---|
[11322] | 1260 | !! |
---|
| 1261 | !!========================================================= |
---|
| 1262 | !! |
---|
| 1263 | !! Hanning filter |
---|
| 1264 | smth_a = 1._wp / 8._wp |
---|
| 1265 | smth_b = 1._wp / 4._wp |
---|
| 1266 | ! |
---|
[14215] | 1267 | DO_2D( 1, 0, 1, 1 ) |
---|
[12353] | 1268 | zumsk = msk(ji,jj) * msk(ji+1,jj) |
---|
| 1269 | zdX ( ji, jj ) = ( pvar2d( ji+1,jj ) - pvar2d( ji ,jj ) ) * zumsk |
---|
| 1270 | END_2D |
---|
[13214] | 1271 | |
---|
[14215] | 1272 | DO_2D( 1, 1, 1, 0 ) |
---|
[12353] | 1273 | zvmsk = msk(ji,jj) * msk(ji,jj+1) |
---|
| 1274 | zdY ( ji, jj ) = ( pvar2d( ji, jj+1 ) - pvar2d( ji ,jj ) ) * zvmsk |
---|
[13214] | 1275 | END_2D |
---|
| 1276 | |
---|
[14215] | 1277 | DO_2D( 0, 0, 1, 0 ) |
---|
[12353] | 1278 | zFY ( ji, jj ) = zdY ( ji, jj ) & |
---|
| 1279 | & + smth_a* ( (zdX ( ji, jj+1 ) - zdX( ji-1, jj+1 )) & |
---|
| 1280 | & - (zdX ( ji, jj ) - zdX( ji-1, jj )) ) |
---|
[13214] | 1281 | END_2D |
---|
[11322] | 1282 | |
---|
[14215] | 1283 | DO_2D( 1, 0, 0, 0 ) |
---|
[12353] | 1284 | zFX( ji, jj ) = zdX( ji, jj ) & |
---|
| 1285 | & + smth_a*( (zdY( ji+1, jj ) - zdY( ji+1, jj-1)) & |
---|
| 1286 | & - (zdY( ji , jj ) - zdY( ji , jj-1)) ) |
---|
| 1287 | END_2D |
---|
[11322] | 1288 | |
---|
[13295] | 1289 | DO_2D( 0, 0, 0, 0 ) |
---|
[12353] | 1290 | pvar2d( ji ,jj ) = pvar2d( ji ,jj ) & |
---|
| 1291 | & + msk(ji,jj) * smth_b * ( & |
---|
| 1292 | & zFX( ji, jj ) - zFX( ji-1, jj ) & |
---|
| 1293 | & +zFY( ji, jj ) - zFY( ji, jj-1 ) ) |
---|
| 1294 | END_2D |
---|
[13214] | 1295 | |
---|
[11322] | 1296 | !--------------------------------------------------------------------------------------------------- |
---|
| 1297 | END SUBROUTINE smooth_pblh |
---|
| 1298 | !=================================================================================================== |
---|
| 1299 | |
---|
[11305] | 1300 | !!====================================================================== |
---|
| 1301 | END MODULE ablmod |
---|