MODULE pbl_parameters_mod !======================================================================================================================! ! Subject: !--------- ! Module used to compute the friction velocity, temperature, and monin_obukhov length from temperature, wind field and thermals outputs. ! Also interpolate theta and u;v in the surface layer based on the Monin Obukhov theory ! These are diagnostics only and do not influence the code. !----------------------------------------------------------------------------------------------------------------------! ! Reference: !----------- ! Colaïtis, A., Spiga, A., Hourdin, F., Rio, C., Forget, F., and Millour, E. (2013), A thermal plume model for the Martian convective boundary layer, J. Geophys. Res. Planets, 118, 1468–1487, doi:10.1002/jgre.20104. ! !======================================================================================================================! IMPLICIT NONE CONTAINS SUBROUTINE pbl_parameters(ngrid,nlay,ps,pplay,pz0,pg,zzlay,zzlev,pu,pv,wstar_in,hfmax,zmax,tke,pts,ph,pqvap,pqsurf,mumean,z_out,n_out, & T_out,u_out,ustar,tstar,vhf,vvv) USE comcstfi_h use write_output_mod, only: write_output use turb_mod, only: turb_resolved use lmdz_atke_turbulence_ini, only : smmin, ric, cinf, cepsilon, pr_slope, pr_asym, pr_neut, ri0, ri1, cn, rpi use watersat_mod, only: watersat use paleoclimate_mod, only: include_waterbuoyancy IMPLICIT NONE !======================================================================= ! ! Anlysis of the PBL from input temperature, wind field and thermals outputs: compute the friction velocity, temperature, and monin_obukhov length ! and interpolate the potential temperature and winds in the surface layer using the Monin Obukhov theory ! ! ------- ! ! Author: Arnaud Colaitis 09/01/12; adapted in F90 by Lucas Lange 12/2023 ! ------- ! ! Arguments: ! ---------- ! ! inputs: ! ------ ! ngrid size of the horizontal grid ! nlay size of the vertical grid ! ps(ngrid) surface pressure (Pa) ! pplay(ngrid,nlay) pressure levels ! pz0(ngrid) surface roughness length (m) ! pg gravity (m s -2) ! zzlay(ngrid,nlay) height of mid-layers (m) ! zzlev(ngrid,nlay+1) height of layers interfaces (m) ! pu(ngrid,nlay) u component of the wind (m/s) ! pv(ngrid,nlay) v component of the wind (m/s) ! wstar_in(ngrid) free convection velocity in PBL (m/s) ! hfmax(ngrid) maximum vertical turbulent heat flux in thermals (W/m^2) ! zmax(ngrid) height reached by the thermals (pbl height) (m) ! tke(ngrid,nlay+1) turbulent kinetic energy (J/kg) ! pts(ngrid) surface temperature (K) ! ph(ngrid,nlay) potential temperature T*(p/ps)^kappa (k) ! z_out(n_out) heights of interpolation (m) ! n_out number of points for interpolation (1) ! ! outputs: ! ------ ! ! T_out(ngrid,n_out) interpolated potential temperature on z_out (K) ! u_out(ngrid,n_out) interpolated winds on z_out (m/s) ! ustar(ngrid) friction velocity (m/s) ! tstar(ngrid) friction temperature (K) ! ! !======================================================================= ! !----------------------------------------------------------------------- ! Declarations: ! ------------- #include "callkeys.h" ! Arguments: ! ---------- INTEGER, INTENT(IN) :: ngrid,nlay,n_out ! size of the horizontal and vertical grid, interpolated altitudes for the surface layer REAL, INTENT(IN) :: pz0(ngrid) ! surface roughness REAL, INTENT(IN) :: ps(ngrid),pplay(ngrid,nlay) ! surface pressure and pressure levels REAL, INTENT(IN) :: pg,zzlay(ngrid,nlay),zzlev(ngrid,nlay) ! surface gravity, altitude of the interface and mid-layers REAL, INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay) ! winds REAL, INTENT(IN) :: wstar_in(ngrid),hfmax(ngrid),zmax(ngrid) ! free convection velocity , vertical turbulent heat, pbl height REAL, INTENT(IN) :: tke(ngrid,nlay+1) ! Turbulent kinetic energy REAL, INTENT(IN) :: pts(ngrid),ph(ngrid,nlay) ! surface temperature, potentiel temperature REAL, INTENT(IN) :: z_out(n_out) ! altitudes of the interpolation in the surface layer REAL, INTENT(IN) :: mumean(ngrid) ! Molecular mass of the atmosphere (kg/mol) REAL, INTENT(IN) :: pqvap(ngrid,nlay) ! Water vapor mixing ratio (kg/kg) to account for vapor flottability REAL, INTENT(IN) :: pqsurf(ngrid) ! Water ice frost on the surface (kg/m^2) to account for vapor flottability ! Outputs: ! -------- REAL, INTENT(OUT) :: T_out(ngrid,n_out),u_out(ngrid,n_out) ! Temperature and wind of the interpolation in the surface layer REAL, INTENT(OUT) :: ustar(ngrid), tstar(ngrid) ! Friction velocity and temperature ! Local: ! ------ INTEGER ig,k,n ! Loop variables INTEGER ii(1) ! Index to compute the pbl mixed layer temperature REAL karman,nu ! Von Karman constant and fluid kinematic viscosity DATA karman,nu/.41,0.001/ !$OMP THREADPRIVATE(karman,nu) SAVE karman,nu REAL zout ! altitude to interpolate (local variable during the loop on z_out) (m) REAL rib(ngrid) ! Bulk Richardson number (1) REAL fm(ngrid) ! stability function for momentum (1) REAL fh(ngrid) ! stability function for heat (1) REAL z1z0,z1z0t ! ratios z1/z0 and z1/z0T (1) ! Formalism for stability functions fm= 1/phim^2; fh = 1/(phim*phih) where ! phim = 1+betam*zeta or (1-bm*zeta)**am ! phih = alphah + betah*zeta or alphah(1.-bh*zeta)**ah ! ah and am are assumed to be -0.25 and -0.5 respectively ! lambda is an intermediate variable to simplify the expression REAL betam, betah, alphah, bm, bh, lambda REAL cdn(ngrid),chn(ngrid) ! neutral momentum and heat drag coefficient (1) REAL pz0t ! initial thermal roughness length z0t for the iterative scheme (m) REAL ric_colaitis ! critical richardson number (1) REAL reynolds(ngrid) ! Reynolds number (1) REAL prandtl(ngrid) ! Prandtl number (1) REAL pz0tcomp(ngrid) ! Computed z0t during the iterative scheme (m) REAL ite ! Number of iteration (1) REAL residual ! Residual between pz0t and pz0tcomp (m) REAL itemax ! Maximum number of iteration (1) REAL tol_iter ! Tolerance for the residual to ensure convergence (1) REAL pcdv(ngrid),pcdh(ngrid) ! momentum and heat drag coefficient (1) REAL zu2(ngrid) ! Near surface winds (m/s) REAL x(ngrid) ! z/zi (1) REAL dvhf(ngrid),dvvv(ngrid) ! dimensionless vertical heat flux and dimensionless vertical velocity variance REAL vhf(ngrid),vvv(ngrid) ! vertical heat flux (W/m^2) and vertical velocity variance (m/s) REAL Teta_out(ngrid,n_out) ! Temporary variable to compute interpolated potential temperature (K) REAL sm ! Stability function computed with the ATKE scheme REAL f_ri_cd_min ! minimum of the stability function with the ATKE scheme REAL ric_4interp ! critical richardson number which is either the one from Colaitis et al. (2013) or ATKE (1) REAL tsurf_v(ngrid) ! Virtual surface temperature, accounting for vapor flottability REAL temp_v(ngrid) ! Potential virtual air temperature in the frist layer, accounting for vapor flottability REAL mu_h2o ! Molecular mass of water vapor (kg/mol) REAL tol_frost ! Tolerance to consider the effect of frost on the surface REAL qsat(ngrid) ! saturated mixing ratio of water vapor ! Code: ! -------- !------------------------------------------------------------------------ !------------------------------------------------------------------------ ! PART I : RICHARDSON/REYNOLDS/THERMAL_ROUGHNESS/STABILITY_COEFFICIENTS !------------------------------------------------------------------------ !------------------------------------------------------------------------ ! Initialisation : ustar(:)=0. tstar(:)=0. reynolds(:)=0. pz0t = 0. pz0tcomp(:) = 0.1*pz0(:) rib(:)=0. pcdv(:)=0. pcdh(:)=0. itemax= 10 tol_iter = 0.01 f_ri_cd_min = 0.01 mu_h2o = 18e-3 tol_frost = 1e-4 tsurf_v(:) = 0. temp_v(:) = 0. ! this formulation assumes alphah=1., implying betah=betam ! We use Dyer et al. parameters, as they cover a broad range of Richardson numbers : bm = 16. bh = 16. alphah = 1. betam = 5. betah = 5. lambda = (sqrt(bh/bm))/alphah ric_colaitis = betah/(betam**2) if(callatke) then ric_4interp = ric else ric_4interp = ric_colaitis endif IF(include_waterbuoyancy) then temp_v(:) = ph(:,1)*(1.+pqvap(:,1)/(mu_h2o/mumean(:)))/(1.+pqvap(:,1)) DO ig = 1,ngrid IF(pqsurf(ig).gt.tol_frost) then call watersat(1,pts(ig),pplay(ig,1),qsat(ig)) tsurf_v(ig) = pts(ig)*(1.+qsat(ig)/(mu_h2o/mumean(ig)))/(1.+qsat(ig)) ELSE tsurf_v(ig) = pts(ig)*(1.+pqvap(ig,1)/(mu_h2o/mumean(ig)))/(1.+pqvap(ig,1)) ENDIF ENDDO ELSE tsurf_v(:) = pts(:) temp_v(:) = ph(:,1) ENDIF DO ig=1,ngrid ite = 0. residual = abs(pz0tcomp(ig)-pz0t) z1z0 = zzlay(ig,1)/pz0(ig) cdn(ig) = karman/log(z1z0) cdn(ig) = cdn(ig)*cdn(ig) zu2(ig) = pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1) + (log(1.+0.7*wstar_in(ig) + 2.3*wstar_in(ig)**2))**2 ! near surface winds, accounting for buyoncy ! IF(turb_resolved) zu2(ig)=MAX(zu2(ig),1.) DO WHILE((residual .gt. tol_iter*pz0(ig)) .and. (ite .lt. itemax)) ! Computations of z0T; iterated until z0T converges pz0t = pz0tcomp(ig) z1z0t=zzlay(ig,1)/pz0t chn(ig) = cdn(ig)*log(z1z0)/log(z1z0t) IF ((pu(ig,1) .ne. 0.) .or. (pv(ig,1) .ne. 0.)) then ! Richardson number formulation proposed by D.E. England et al. (1995) rib(ig) = (pg/tsurf_v(ig))*sqrt(zzlay(ig,1)*pz0(ig))*(((log(zzlay(ig,1)/pz0(ig)))**2)/(log(zzlay(ig,1)/pz0t)))*(temp_v(ig)-tsurf_v(ig))/zu2(ig) ELSE print*,'warning, infinite Richardson at surface' print*,pu(ig,1),pv(ig,1) rib(ig) = ric_colaitis ENDIF ! Compute the stability functions fm; fh depending on the stability of the surface layer IF(callatke) THEN ! Computation following parameterizaiton from ATKE IF(rib(ig) .gt. 0) THEN ! STABLE BOUNDARY LAYER : sm = MAX(smmin,cn*(1.-rib(ig)/ric)) ! prandlt expression from venayagamoorthy and stretch 2010, Li et al 2019 prandtl(ig) = pr_neut*exp(-pr_slope/pr_neut*rib(ig)+rib(ig)/pr_neut) + rib(ig) * pr_slope fm(ig) = max((sm**(3./2.) * sqrt(cepsilon) * (1 - rib(ig) / prandtl(ig))**(1./2.)),f_ri_cd_min) fh(ig) = max((fm(ig) / prandtl(ig)), f_ri_cd_min) ELSE ! UNSTABLE BOUNDARY LAYER : sm = 2./rpi * (cinf - cn) * atan(-rib(ig)/ri0) + cn prandtl(ig) = -2./rpi * (pr_asym - pr_neut) * atan(rib(ig)/ri1) + pr_neut fm(ig) = MAX((sm**(3./2.) * sqrt(cepsilon) * (1 - rib(ig) / prandtl(ig))**(1./2.)),f_ri_cd_min) fh(ig) = MAX((fm(ig) / prandtl(ig)), f_ri_cd_min) ENDIF ! Rib < 0 ELSE ! Computation following parameterizaiton from from D.E. England et al. (95) IF (rib(ig) .gt. 0.) THEN ! STABLE BOUNDARY LAYER : prandtl(ig) = 1. IF(rib(ig) .lt. ric_colaitis) THEN ! Assuming alphah=1. and bh=bm for stable conditions : fm(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2 fh(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2 ELSE ! For Ri>Ric, we consider Ri->Infinity => no turbulent mixing at surface fm(ig) = 1. fh(ig) = 1. ENDIF ELSE ! UNSTABLE BOUNDARY LAYER : fm(ig) = sqrt(1.-lambda*bm*rib(ig)) fh(ig) = (1./alphah)*((1.-lambda*bh*rib(ig))**0.5)*(1.-lambda*bm*rib(ig))**0.25 prandtl(ig) = alphah*((1.-lambda*bm*rib(ig))**0.25)/((1.-lambda*bh*rib(ig))**0.5) ENDIF ! Rib < 0 ENDIF ! atke ! Recompute the reynolds and z0t reynolds(ig) = karman*sqrt(fm(ig))*sqrt(zu2(ig))*pz0(ig)/(log(z1z0)*nu) pz0tcomp(ig) = pz0(ig)*exp(-karman*7.3*(reynolds(ig)**0.25)*(prandtl(ig)**0.5)+5*karman) residual = abs(pz0t-pz0tcomp(ig)) ite = ite+1 ENDDO ! of while ! Compute the coefficients Cdv; Cdh : neutral coefficient x stability functions pcdv(ig) = cdn(ig)*fm(ig) pcdh(ig) = chn(ig)*fh(ig) pz0t = 0. ! for next grid point ENDDO ! of ngrid !------------------------------------------------------------------------ !------------------------------------------------------------------------ ! PART II : USTAR/TSTAR/U_OUT/TETA_OUT COMPUTATIONS !------------------------------------------------------------------------ !------------------------------------------------------------------------ ! u* theta* computation DO n=1,n_out zout = z_out(n) DO ig=1,ngrid IF (rib(ig) .ge. ric_4interp) THEN !stable case, no turbulence ustar(ig) = 0. tstar(ig) = 0. ELSE ustar(ig) = sqrt(pcdv(ig))*sqrt(zu2(ig)) ! By definition, u* = sqrt(tau/U^2) = sqrt(Cdv) tstar(ig) = -pcdh(ig)*(pts(ig)-ph(ig,1))/sqrt(pcdv(ig)) ! theta* = (T1-Tsurf)*Cdh(ig)/u* ENDIF ! Interpolation: IF(zout .lt. pz0tcomp(ig)) THEN ! If z_out is in the surface layer , theta = tsurf; u = 0 u_out(ig,n) = 0. Teta_out(ig,n) = pts(ig) ELSEIF (zout .lt. pz0(ig)) THEN u_out(ig,n) = 0. ELSE IF (rib(ig) .ge. ric_4interp) THEN ! ustar=tstar=0 (and fm=fh=0) because no turbulence u_out(ig,n) = 0 Teta_out(ig,n) = pts(ig) ELSE ! We use the MO theory: du/dz = u*/(kappa z) *1 /sqrt(fm); dtheta/zd = theta*/(Pr kappa z) * fm/fh u_out(ig,n)= ustar(ig)*log(zout/pz0(ig))/(karman*sqrt(fm(ig))) Teta_out(ig,n)=pts(ig)+(tstar(ig)*sqrt(fm(ig))*log(zout/(pz0tcomp(ig)))/(karman*fh(ig))) ENDIF ENDIF ENDDO ! loop on n_out ! when using convective adjustment without thermals, a vertical potential temperature ! profile is assumed up to the thermal roughness length. Hence, theoretically, theta ! interpolated at any height in the surface layer is theta at the first level. IF ((.not.calltherm).and.(calladj)) THEN Teta_out(:,n)=ph(:,1) u_out(:,n)=(sqrt(cdn(:))*sqrt(zu2(ig))/karman)*log(zout/pz0(:)) ENDIF T_out(:,n) = Teta_out(:,n)*(exp((zout/zzlay(:,1))*(log(pplay(:,1)/ps))))**rcp ! not sure of what is done here: hydrostatic correction to account for the difference of pressure between zout and z1 ? ENDDO !of n=1,n_out !------------------------------------------------------------------------ !------------------------------------------------------------------------ ! PART III : VERTICAL_VELOCITY_VARIANCE/VERTICAL_TURBULENT_FLUX PROFILES !------------------------------------------------------------------------ !------------------------------------------------------------------------ ! We follow Spiga et. al 2010 (QJRMS) ! ------------ IF (calltherm) THEN DO ig=1, ngrid IF (zmax(ig) .gt. 0.) THEN x(ig) = zout/zmax(ig) ELSE x(ig) = 999. ENDIF ENDDO DO ig=1, ngrid ! dimensionless vertical heat flux IF (x(ig) .le. 0.3) THEN dvhf(ig) = ((-3.85/log(x(ig)))+0.07*log(x(ig))) *exp(-4.61*x(ig)) ELSEIF (x(ig) .le. 1.) THEN dvhf(ig) = -1.52*x(ig) + 1.24 ELSE dvhf(ig) = 0. ENDIF ! dimensionless vertical velocity variance IF (x(ig) .le. 1.) THEN dvvv(ig) = 2.05*(x(ig)**(2./3.))*(1.-0.64*x(ig))**2 ELSE dvvv(ig) = 0. ENDIF ENDDO vhf(:) = dvhf(:)*hfmax(:) vvv(:) = dvvv(:)*(wstar_in(:))**2 ENDIF ! of if calltherm !------------------------------------------------------------------------ !------------------------------------------------------------------------ ! PART IV : Outputs !------------------------------------------------------------------------ !------------------------------------------------------------------------ #ifndef MESOSCALE call write_output('tke_pbl','Tke in the pbl after physic','J/kg',tke(:,:)) call write_output('rib_pbl','Richardson in pbl parameter','m/s',rib(:)) call write_output('cdn_pbl','neutral momentum coef','m/s',cdn(:)) call write_output('fm_pbl','momentum stab function','m/s',fm(:)) call write_output('uv','wind norm first layer','m/s',sqrt(zu2(:))) call write_output('uvtrue','wind norm first layer','m/s',sqrt(pu(:,1)**2+pv(:,1)**2)) call write_output('chn_pbl','neutral momentum coef','m/s',chn(:)) call write_output('fh_pbl','momentum stab function','m/s',fh(:)) call write_output('B_pbl','buoyancy','m/',(ph(:,1)-pts(:))/pts(:)) call write_output('zot_pbl','buoyancy','ms',pz0tcomp(:)) call write_output('zz1','1st layer altitude','m',zzlay(:,1)) #endif RETURN END SUBROUTINE pbl_parameters END MODULE pbl_parameters_mod