MODULE vdif_cd_mod !======================================================================================================================! ! Subject: !--------- ! Module used to compute the exchange coefficient in the surface layers Cd; Ch !----------------------------------------------------------------------------------------------------------------------! ! 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 vdif_cd(ngrid,nlay,nslope,pz0,pg,pz,pp,pu,pv,wstar,pts,ph,mumean,pqvap,pqsurf,write_outputs,pcdv,pcdh) use comcstfi_h use turb_mod, only: turb_resolved use watersat_mod, only: watersat use lmdz_atke_turbulence_ini, only : smmin, ric, cinf, cepsilon, pr_slope, pr_asym, pr_neut, ri0, ri1, cn, rpi use paleoclimate_mod, only: include_waterbuoyancy use write_output_mod, only: write_output use comslope_mod, ONLY: iflat IMPLICIT NONE include "callkeys.h" !======================================================================= ! ! Subject: computation of the surface drag coefficient using the ! ------- approch developed by Loui for ECMWF. ! ! Author: Frederic Hourdin 15 /10 /93 ! Modified by : Arnaud Colaitis 03/08/11; rewritten by LL to switch to F90 ! ------- ! ! Arguments: ! ---------- ! ! inputs: ! ------ ! ngrid size of the horizontal grid ! pg gravity (m s -2) ! pz(ngrid,nlay) height of layers ! pp(ngrid,nlay) pressure of layers ! pu(ngrid,nlay) u component of the wind ! pv(ngrid,nlay) v component of the wind ! pts(ngrid) surface temperature ! ph(ngrid) potential temperature T*(p/ps)^kappa ! mumean Molecular mass of the atmosphere (kg/mol) ! pqvap(ngrid,nlay) Water vapor mixing ratio (kg/kg) to account for vapor flottability ! pqsurf(ngrid) Water ice frost on the surface(kg/m^2) to account for vapor flottability ! ! outputs: ! -------- ! pcdv(ngrid) Cd for the wind ! pcdh(ngrid) Cd for potential temperature ! !======================================================================= !----------------------------------------------------------------------- ! Declarations: ! ------------- ! Arguments: ! Inputs: ! ---------- INTEGER, INTENT(IN) :: ngrid,nlay,nslope ! Number of points in the physical grid and atmospheric grid REAL, INTENT(IN) :: pz0(ngrid) ! Surface Roughness (m) REAL, INTENT(IN) :: pg ! Mars gravity (m/s^2) REAL, INTENT(IN) :: pz(ngrid,nlay),pp(ngrid,nlay) ! Layers pseudo altitudes (m) and pressure (Pa) REAL, INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay) ! Zonal and Meriditionnal winds (m/s) REAL, INTENT(IN) :: pts(ngrid,nslope),ph(ngrid,nlay) ! Surface Temperature and atmospheric temperature (K) REAL, INTENT(IN) :: wstar(ngrid) ! Vertical velocity due to turbulence (m/s) 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,nslope) ! Water ice frost on the surface (kg/m^2) to account for vapor flottability LOGICAL, INTENT(IN) :: write_outputs ! write_output in xios or not. ! Outputs: ! ---------- REAL, INTENT(OUT) :: pcdv(ngrid,nslope) ! momentum drag coefficient (1) REAL, INTENT(OUT) :: pcdh(ngrid,nslope) ! heat drag coefficient (1) ! Local: ! ------ INTEGER ig,islope ! Loop variable REAL karman,nu ! Von Karman constant and fluid kinematic viscosity LOGICAL firstcal DATA karman,nu/.41,0.001/ DATA firstcal/.true./ SAVE karman,nu !$OMP THREADPRIVATE(karman,nu) REAL rib(ngrid,nslope) ! Bulk Richardson number [virtual or dry] (1) REAL rib_dry(ngrid,nslope) ! Dry bulk Richardson number [virtual or dry] (1) REAL fm(ngrid,nslope) ! stability function for momentum (1) REAL fh(ngrid,nslope) ! stability function for heat (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 pz0t ! initial thermal roughness length z0t for the iterative scheme (m) REAL ric_colaitis ! critical richardson number (1) REAL reynolds(ngrid,nslope) ! Reynolds number (1) REAL prandtl(ngrid) ! Prandtl number (1) REAL pz0tcomp(ngrid) ! computed z0t during the iterative scheme (m) REAL z0t(ngrid,nslope) ! computed z0t at the last step the iterative scheme (m) REAL ite ! Number of iteration (1) REAL itemax ! Maximum number of iteration (1) REAL residual ! Residual between pz0t and pz0tcomp (m) REAL tol_iter ! Tolerance for the residual to ensure convergence (1= REAL zu2(ngrid) ! Near surface winds (m/s) REAL cdn(ngrid) ! neutral momentum drag coefficient (1) REAL chn(ngrid) ! neutral heat drag coefficient (1) REAL z1z0,z1z0t ! ratios z1/z0 and z1/z0T REAL z1,zcd0 ! Neutral roughness (m) and Cd/Ch coefficient when call richls is not called REAL tsurf_v(ngrid,nslope) ! 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 REAL sm ! Stability function computed with the ATKE scheme REAL f_ri_cd_min ! minimum of the stability function with the ATKE scheme ! Code: ! -------- ! Initialisation itemax= 10 tol_iter = 0.01 mu_h2o = 18e-3 tol_frost = 1e-4 reynolds(:,:) = 0. pz0t = 0. pz0tcomp(:) = 0.1*pz0(:) rib(:,:) = 0. rib_dry(:,:) = 0. pcdv(:,:) = 0. pcdh(:,:) = 0. z0t(:,:) = 0. fm(:,:) = 0. fh(:,:) = 0. f_ri_cd_min = 0.01 ! 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(include_waterbuoyancy) then temp_v(:) = ph(:,1)*(1.+pqvap(:,1)/(mu_h2o/mumean(:)))/(1.+pqvap(:,1)) DO islope = 1,nslope DO ig = 1,ngrid IF(pqsurf(ig,islope).gt.tol_frost) then call watersat(1,pts(ig,islope),pp(ig,1),qsat(ig)) tsurf_v(ig,islope) = pts(ig,islope)*(1.+qsat(ig)/(mu_h2o/mumean(ig)))/(1.+qsat(ig)) ELSE tsurf_v(ig,islope) = pts(ig,islope)*(1.+pqvap(ig,1)/(mu_h2o/mumean(ig)))/(1.+pqvap(ig,1)) ENDIF ENDDO ENDDO ELSE tsurf_v(:,:) = pts(:,:) temp_v(:) = ph(:,1) ENDIF IF(.not.callrichsl) then ! Original formulation as in LMDZ Earth: Cd = Ch = (kappa/(ln(1+z1/z0)))^2 ! NB: In forget et al., 1999, it's Cd = Ch = (kappa/(ln(z1/z0)))^2 DO ig = 1,ngrid z1 = 1.E+0 + pz(ig,1)/pz0(ig) zcd0 = karman/log(z1) zcd0 = zcd0*zcd0 DO islope = 1,nslope pcdv(ig,islope) = zcd0 pcdh(ig,islope) = zcd0 ENDDO ENDDO ELSE ! We follow the parameterization from Colaitis et al., 2013 (supplementary material) DO islope = 1,nslope DO ig=1,ngrid ite = 0. residual = 100*tol_iter*pz0(ig) z1z0=pz(ig,1)/pz0(ig) cdn(ig)=karman/log(z1z0) cdn(ig)=cdn(ig)*cdn(ig) DO WHILE((residual .gt. tol_iter*pz0(ig)) .and. (ite .lt. itemax)) ! Computations of z0T; iterated until z0T converges pz0t = pz0tcomp(ig) z1z0t=pz(ig,1)/pz0t chn(ig)=cdn(ig)*log(z1z0)/log(z1z0t) IF ((pu(ig,1) .ne. 0.) .or. (pv(ig,1) .ne. 0.)) then zu2(ig) = pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1) + (log(1.+0.7*wstar(ig) + 2.3*wstar(ig)**2))**2 ! Near surface winds account for buoyancy IF(turb_resolved) zu2(ig)=MAX(zu2(ig),1.) ! Richardson number formulation proposed by D.E. England et al. (1995) rib(ig,islope) = (pg/tsurf_v(ig,islope))*sqrt(pz(ig,1)*pz0(ig))*(((log(pz(ig,1)/pz0(ig)))**2)/(log(pz(ig,1)/pz0t)))*(temp_v(ig)-tsurf_v(ig,islope))/zu2(ig) rib_dry(ig,islope) = (pg/pts(ig,islope))*sqrt(pz(ig,1)*pz0(ig))*(((log(pz(ig,1)/pz0(ig)))**2)/(log(pz(ig,1)/pz0t)))*(ph(ig,1)-pts(ig,islope))/zu2(ig) ELSE print*,'warning, infinite Richardson at surface' print*,pu(ig,1),pv(ig,1) rib(ig,islope) = ric_colaitis rib_dry(ig,islope) = 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,islope) .gt. 0) THEN ! STABLE BOUNDARY LAYER :include_waterbuoyancy sm = MAX(smmin,cn*(1.-rib(ig,islope)/ric)) ! prandlt expression from venayagamoorthy and stretch 2010, Li et al 2019 prandtl(ig) = pr_neut*exp(-pr_slope/pr_neut*rib(ig,islope)+rib(ig,islope)/pr_neut) + rib(ig,islope) * pr_slope fm(ig,islope) = max((sm**(3./2.) * sqrt(cepsilon) * (1 - rib(ig,islope) / prandtl(ig))**(1./2.)),f_ri_cd_min) fh(ig,islope) = max((fm(ig,islope) / prandtl(ig)), f_ri_cd_min) ELSE ! UNSTABLE BOUNDARY LAYER : sm = 2./rpi * (cinf - cn) * atan(-rib(ig,islope)/ri0) + cn prandtl(ig) = -2./rpi * (pr_asym - pr_neut) * atan(rib(ig,islope)/ri1) + pr_neut fm(ig,islope) = MAX((sm**(3./2.) * sqrt(cepsilon) * (1 - rib(ig,islope) / prandtl(ig))**(1./2.)),f_ri_cd_min) fh(ig,islope) = MAX((fm(ig,islope) / prandtl(ig)), f_ri_cd_min) ENDIF ! Rib < 0 ELSE ! Computation following parameterizaiton from from D.E. England et al. (95) IF (rib(ig,islope) .gt. 0.) THEN ! STABLE BOUNDARY LAYER : prandtl(ig) = 1. IF(rib(ig,islope) .lt. ric_colaitis) THEN ! Assuming alphah=1. and bh=bm for stable conditions : fm(ig,islope) = ((ric_colaitis-rib(ig,islope))/ric_colaitis)**2 fh(ig,islope) = ((ric_colaitis-rib(ig,islope))/ric_colaitis)**2 ELSE ! For Ri>Ric, we consider Ri->Infinity => no turbulent mixing at surface fm(ig,islope) = 1. fh(ig,islope) = 1. ENDIF ELSE ! UNSTABLE BOUNDARY LAYER : fm(ig,islope) = sqrt(1.-lambda*bm*rib(ig,islope)) fh(ig,islope) = (1./alphah)*((1.-lambda*bh*rib(ig,islope))**0.5)*(1.-lambda*bm*rib(ig,islope))**0.25 prandtl(ig) = alphah*((1.-lambda*bm*rib(ig,islope))**0.25)/((1.-lambda*bh*rib(ig,islope))**0.5) ENDIF ! Rib < 0 ENDIF ! atke ! Recompute the reynolds and z0t reynolds(ig,islope) = karman*sqrt(fm(ig,islope))*sqrt(zu2(ig))*pz0(ig)/(log(z1z0)*nu) pz0tcomp(ig) = pz0(ig)*exp(-karman*7.3*(reynolds(ig,islope)**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,islope)=cdn(ig)*fm(ig,islope) pcdh(ig,islope)=chn(ig)*fh(ig,islope) z0t(ig,islope) = pz0tcomp(ig) pz0t = 0. ! for next grid point ENDDO ! of ngrid enddo ENDIF !of if call richardson surface layer IF(write_outputs) then call write_output("rib_dry_vdif_cd","Dry Richardson number in vdif_cd","1",rib_dry(:,iflat)) call write_output("rib_vdif_cd","Richardson number in vdif_cd","1",rib(:,iflat)) call write_output("fm_vdif_cd","Momentum stability function in vdif_cd","1",fm(:,iflat)) call write_output("fh_vdif_cd","Heat stability function in vdif_cd","1",fh(:,iflat)) call write_output("z0t_vdif_cd","Thermal roughness length in vdif_cd","m",z0t(:,iflat)) call write_output("z0_vdif_cd","Momentum roughness length in vdif_cd","m",pz0(:)) call write_output("Reynolds_vdif_cd","Reynolds number in vdif_cd","1",reynolds(:,iflat)) ENDIF RETURN END SUBROUTINE vdif_cd END MODULE vdif_cd_mod