module kphysmod 5
implicit none
contains
subroutine radiation(topflux,acoef,nfield) 1,2
use kvar
; use kcnst
real(mrl) :: colwater,topflux,acoef
integer nfield
radflux(nf)=topflux
colwater=0.
if (nfield.eq.1) then
do k=nf-1,1,-1
colwater=colwater+(mf(k+1)-mf(k))*(qc(k)+qr(k))
radflux(k)=radflux(nf)*exp(-colwater*acoef)
radrate(k)=-(radflux(k+1)-radflux(k))/(mf(k+1)-mf(k))/1004.
end do
else if (nfield.eq.2) then
do k=nf-1,1,-1
colwater=colwater+(mf(k+1)-mf(k))*smoke(k)
radflux(k)=radflux(nf)*exp(-colwater*acoef)
radrate(k)=-(radflux(k+1)-radflux(k))/(mf(k+1)-mf(k))/1004.
end do
end if
end subroutine radiation
subroutine blgrid_prep 3,4
use kvar
; use kcnst
; use blgridmod
real(mrl) hite1,hite2,hite3,hite4
!real(mrl) blgsq,blgaq,blgalf,blgbeta
real(mrl),allocatable :: zb(:)
integer nit
!blgsq=.2; blgalf=2; blgbeta=2
allocate(zb(nf))
hite1=150.*.5
hite2=hite-300*.5
hite3=hite+100*.5
hite4=hite3+100.*.5+.4*max(hite3-zf(1),0._mrl)
nit=4
if (time.eq.0.) then
nit=60
blgaq=blgsq
end if
do n=1,nit
call blgrid
(zfdef,zb,zfnext,hite1,hite2,hite3,hite4,&
blgsq,blgaq,blgalf,blgbeta)
end do
deallocate(zb)
end subroutine blgrid_prep
subroutine nocturnal_grid 2,2
use kvar
; use kcnst
real(mrl) dum,shrink
integer nexp
zf(1)=zbottom
zf(nf)=ztop
nexp=4
shrink=.05
do n=2,nf-1
dum=(n-1.)/(nf-1.)
zf(n) = zf(1) + (zf(nf)-zf(1))*&
(shrink*dum+(2-shrink)*dum**(nexp+1))/(1+dum**nexp)
! write (*,*) dum,zf(n),n
end do
print *, ' nocturnal_grid initialized'
return
end subroutine nocturnal_grid
subroutine uniform_grid 1,2
use kvar
; use kcnst
do i=1,nf
zf(i)=(i-1)*ztop/(nf-1)
end do
print *, ' uniform_grid initialized'
return
end subroutine uniform_grid
subroutine z_coords 5,2
use kvar
; use kcnst
do i=1,nu
zu(i)=0.5*(zf(i)+zf(i+1))
dzf(i)=zf(i+1)-zf(i)
end do
do i=2,nf-1
dzu(i)=zu(i)-zu(i-1)
end do
dzu(1)=zu(1)
dzu(nf)=zf(nf)-zu(nu)
rdzf=1/dzf
rdzu=1/dzu
return
end subroutine z_coords
subroutine init_fields 1,2
use kvar
; use kcnst
ug=ugeo
vg=vgeo
smoke=0.
qr=0.
qc=0.
nc=0.
nr=0.
ccn=0.
smoke_sed=0.
qr_sed=-.05
nr_sed=0.
do n=1,nu
if (zu(n).lt.zinit2) then
th(n)=tinit1+(zu(n)-zinit1)*(tinit2-tinit1)/(zinit2-zinit1)
qv(n)=qinit1+(zu(n)-zinit1)*(qinit2-qinit1)/(zinit2-zinit1)
u(n)=uinit1+(zu(n)-zinit1)*(uinit2-uinit1)/(zinit2-zinit1)
smoke(n)=1.
else if (zu(n).lt.zinit3) then
th(n)=tinit2+(zu(n)-zinit2)*(tinit3-tinit2)/(zinit3-zinit2)
qv(n)=qinit2+(zu(n)-zinit2)*(qinit3-qinit2)/(zinit3-zinit2)
u(n)=uinit2+(zu(n)-zinit2)*(uinit3-uinit2)/(zinit3-zinit2)
smoke(n)=1.-(zu(n)-zinit2)/(zinit3-zinit2)
else if (zu(n).lt.zinit4) then
th(n)=tinit3+(zu(n)-zinit3)*(tinit4-tinit3)/(zinit4-zinit3)
qv(n)=qinit3+(zu(n)-zinit3)*(qinit4-qinit3)/(zinit4-zinit3)
u(n)=uinit3+(zu(n)-zinit3)*(uinit4-uinit3)/(zinit4-zinit3)
else if (zu(n).lt.zinit5) then
th(n)=tinit4+(zu(n)-zinit4)*(tinit5-tinit4)/(zinit5-zinit4)
qv(n)=qinit4+(zu(n)-zinit4)*(qinit5-qinit4)/(zinit5-zinit4)
u(n)=uinit4+(zu(n)-zinit4)*(uinit5-uinit4)/(zinit5-zinit4)
else
th(n)=tinit5+(zu(n)-zinit5)*lapse
qv(n)=qinit5*exp(-.001*(zu(n)-zinit5))
u(n)=uinit5
end if
end do
kzero=0.
if (tkeon.gt.0) tke=0.0003
n=1
do while (zu(n).lt.zinit2)
! smoke(n)=1.
tke(n+1)=tkeinit !new
ccn(n)=50.
n=n+1
end do
!do n=nu,1,-1
! write(*,*) zu(n),th(n),qv(n),ccn(n),smoke(n)
!end do
lengthmix=100.
include '../inc/surface.inc'
print *, ' init_fields successful'
return
end subroutine init_fields
subroutine init_fields_orig ,2
use kvar
; use kcnst
do n=1,nu
! th(n)=280.+min(zu(n),3000._mrl)*lapse
if (zu(n).lt.3000._mrl) then
th(n)=290.+zu(n)*lapse
else
th(n)=290.+3000.*lapse+(zu(n)-3000.)*3.e-3
end if
end do
!th=280.+zu*lapse
!u=1.+zu*lapse
!v=1.+zu*lapse
u=0
v=0
kzero=0.
if (tkeon.gt.0) tke=0.0003
smoke=0.
qr=0.
qv=0.
qc=0.
nc=0.
nr=0.
ccn=0.
smoke_sed=0.
qr_sed=-.05
nr_sed=0.
n=1
th(1)=280.
do while (zu(n).lt.hite)
qv(n)=.003*1. !change 1 to 0 to init as dry
smoke(n)=1.
! qr(n)=qr(1)*.5*(1+cos(6.1415*zu(n)/hite))
ccn(n)=50.
th(n)=th(1)
n=n+1
end do
!do n=1,nu
! write(*,*) th(n)
!end do
!th=th(2,1)
lengthmix=100.
!deallocate(tm)
include '../inc/surface.inc'
print *, ' init_fields successful'
return
end subroutine init_fields_orig
subroutine subsidence 1,3
use kvar
; use kcnst
implicit none
do k=2,nu
massvelf(k)=(mf(k+1)-mf(k-1))/(zf(k+1)-zf(k-1))*f_wbar
(zf(k))
end do
dmfdt(1)=0.
dmfdt(nf)=0.
do k=1,nu
massvelu(k)=0.5*(massvelf(k)+massvelf(k+1))
end do
do k=2,nu
divf(k)=massvelu(k-1)-massvelu(k)
end do
do k=1,nu
divu(k)=massvelf(k)-massvelf(k+1)
end do
do k=1,nu
massvelu(k)=-dmudt(k)+massvelu(k)
end do
do k=2,nu
massvelf(k)=-dmfdt(k)+massvelf(k)
end do
end subroutine subsidence
function f_wbar(zval) 1,2
use kvar
; use kcnst
implicit none
real(mrl) f_wbar,zval
f_wbar=0.
if (zval.lt.wbarhite) then
f_wbar=cos(wbarfreq*2*pi*time)*zval*wbarmag/wbarhite
else
f_wbar=cos(wbarfreq*2*pi*time)*(zf(nf)-zval)&
*wbarmag/(zf(nf)-wbarhite)
end if
return
end function f_wbar
subroutine ad_flux(fq,q,q_sed) 14,2
use kvar
; use kcnst
implicit none
real(mrl),dimension(:),intent(inout) :: fq
real(mrl),dimension(:),intent(in) :: q
real(mrl),dimension(:),intent(in), optional :: q_sed
real(mrl) :: mv
do k=2,nu
mv=massvelf(k)
if (present(q_sed)) then
mv=mv+.5*((mf(k+1)-mf(k))/(zf(k+1)-zf(k))*q_sed(k) &
+(mf(k)-mf(k-1))/(zf(k)-zf(k-1))*q_sed(k-1))
end if
if (mv.lt.0.) then
if (advord.eq.3.and.k.lt.nu) then
fq(k)=fq(k)+mv*(2*q(k-1)+5*q(k)-q(k+1))/6.
else if (advord.eq.2) then
fq(k)=fq(k)+mv*(q(k-1)+q(k))/2.
else
fq(k)=fq(k)+mv*q(k)
end if
else
if (advord.eq.3.and.k.gt.2) then
fq(k)=fq(k)+mv*(2*q(k)+5*q(k-1)-q(k-2))/6.
else if (advord.eq.2) then
fq(k)=fq(k)+mv*(q(k-1)+q(k))/2.
else
fq(k)=fq(k)+mv*q(k-1)
end if
end if
end do
if (present(q_sed)) then
k=1
mv=mv+(mf(k+1)-mf(k))/(zf(k+1)-zf(k))*q_sed(k)
fq(k)=fq(k)+mv*q(k)
end if
end subroutine ad_flux
subroutine ad_flux_f(fq,q) 1,2
use kvar
; use kcnst
implicit none
real(mrl),dimension(:),intent(inout) :: fq
real(mrl),dimension(:),intent(in) :: q
do k=2,nu
if (massvelu(k).lt.0.) then
fq(k)=fq(k)+massvelu(k)*q(k+1)
else
fq(k)=fq(k)+massvelu(k)*q(k)
end if
end do
end subroutine ad_flux_f
subroutine tendencies(dqdt,fq,q) 13,2
use kvar
; use kcnst
implicit none
real(mrl),dimension(:),intent(in) :: fq,q
real(mrl),dimension(:),intent(inout) :: dqdt
do k=1,nu
dqdt(k)=rdmf(k)*(fq(k)-fq(k+1)-divu(k)*q(k))+dqdt(k)
end do
end subroutine tendencies
subroutine mono_corr(fq,q) 8,3
use kvar
; use kcnst
implicit none
real(mrl),dimension(:),intent(inout) :: fq,q
real(mrl) pcheck,qmax,qmin,qk
if (2*(istep/2).lt.istep+10) then
do k=1,nu-1,1
if (k.eq.1) then
qmax=max(q(1),q(2))
qmin=min(q(1),q(2))
else
qmax=max(q(k-1),q(k),q(k+1))
qmin=min(q(k-1),q(k),q(k+1))
end if
qk=dmfold(k)*rdmf(k)*q(k)
pcheck=qk+dt*rdmf(k)*(fq(k)-fq(k+1)-divu(k)*q(k))
if (pcheck.gt.qmax)&
fq(k+1)=-(qmax-qk)/(dt*rdmf(k)) &
+fq(k)-divu(k)*q(k)
if (pcheck.lt.qmin)&
fq(k+1)=-(qmin-qk)/(dt*rdmf(k)) &
+fq(k)-divu(k)*q(k)
end do
else
do k=nu,2,-1
if (k.eq.1) then
qmax=max(q(1),q(2))
qmin=min(q(1),q(2))
else
qmax=max(q(k-1),q(k),q(k+1))
qmin=min(q(k-1),q(k),q(k+1))
end if
pcheck=q(k)+dt*rdmf(k)*(fq(k)-fq(k+1)-divu(k)*q(k))
if (pcheck.gt.qmax)&
fq(k)=(qmax-q(k))/(dt*rdmf(k)) &
+fq(k+1)+divu(k)*q(k)
if (pcheck.lt.qmin)&
fq(k)=(qmin-q(k))/(dt*rdmf(k)) &
+fq(k+1)+divu(k)*q(k)
end do
end if
call positive_corr
(fq,q) !more stringent than above
end subroutine mono_corr
subroutine positive_corr(fq,q) 1,2
use kvar
; use kcnst
implicit none
real(mrl),dimension(:),intent(inout) :: fq,q
real(mrl) pcheck,qmax
do n=2,nu
if (fq(n).gt.0.) then
k=n-1
pcheck=q(k)-dt*rdmf(k)*fq(n)
if (pcheck.lt.0) fq(n)=q(k)/(dt*rdmf(k))
else
k=n
pcheck=q(k)+dt*rdmf(k)*fq(n)
if (pcheck.lt.0) fq(n)=-q(k)/(dt*rdmf(k))
end if
end do
end subroutine positive_corr
subroutine tendencies_f(dqdt,fq,q) 1,2
use kvar
; use kcnst
implicit none
real(mrl),dimension(:),intent(in) :: fq,q
real(mrl),dimension(:),intent(inout) :: dqdt
real(mrl) dummy
do k=2,nu
dqdt(k)=rdmu(k)*(fq(k-1)-fq(k)-divf(k)*q(k))+dqdt(k)
end do
end subroutine tendencies_f
subroutine exner_to_z 1,2
use kvar
; use kcnst
zf(1)=0.
do i=2,nf
zf(i)=zf(i-1)-(cp*th(i-1)/g)*(exner(i)-exner(i-1))
end do
end subroutine exner_to_z
subroutine mass_to_pres_to_exner 1,2
use kvar
; use kcnst
pres=psur-g*mf
exner=(pres/psur)**(rdgas/cp)
end subroutine mass_to_pres_to_exner
subroutine z_to_exner_to_pres 1,2
use kvar
; use kcnst
exner(1)=1.
do i=2,nf
!do i=2,nu
exner(i)=exner(i-1)-g*(zf(i)-zf(i-1))/(cp*th(i-1))
end do
pres=psur*exner**(cp/rdgas)
end subroutine z_to_exner_to_pres
subroutine dzfdt_to_dmfdt 1,2
use kvar
; use kcnst
do k=2,nu
dmfdt(k)=(mf(k+1)-mf(k-1))/(zf(k+1)-zf(k-1))*dzfdt(k)
end do
dmfdt(1)=0.
dmfdt(nf)=0.
do k=1,nu
dmudt(k)=0.5*(dmfdt(k)+dmfdt(k+1))
end do
end subroutine dzfdt_to_dmfdt
subroutine calc_mean_rho 2,2
use kvar
; use kcnst
do i=1,nu
rho(i)=(pres(i)-pres(i+1))*rdzf(i)/g
end do
end subroutine calc_mean_rho
subroutine pres_to_mass 1,2
use kvar
; use kcnst
mf=(psur-pres)/g
end subroutine pres_to_mass
subroutine mass_coords 3,2
use kvar
; use kcnst
dmfold=dmf !temp
do i=1,nu
mu(i)=0.5*(mf(i)+mf(i+1))
dmf(i)=mf(i+1)-mf(i)
end do
do i=2,nf-1
dmu(i)=mu(i)-mu(i-1)
end do
dmu(1)=mu(1)
dmu(nf)=mf(nf)-mu(nu)
rdmu=1/dmu
rdmf=1/dmf
return
end subroutine mass_coords
subroutine diffuse_and_stepx(q2,q1,fq,dqdt,kq) ,3
use kvar
; use kcnst
implicit none
real(mrl), dimension(:), intent(inout) :: q2,q1,fq,dqdt,kq
real(mrl), allocatable, target, dimension(:) :: a,b,c,r
allocate(a(nu),b(nu),c(nu),r(nu))
do i=2,nu
a(i)=-kq(i)*dt*rdzu(i)*rdmf(i)*0.
end do
do i=1,nu-1
c(i)=-kq(i+1)*dt*rdzu(i+1)*rdmf(i)*0.
end do
c(nu)=0.
a(1)=0.
b=1-a-c
!r=q1*dmfold*rdmf+dt*dqdt
r=q1+dt*dqdt
call impl
!call expl
deallocate(a,b,c,r)
contains
subroutine impl 2,2
call tridag
(a,b,c,r,q2)
fq(2:nu)=kq(2:nu)*rdzu(2:nu)*(q2(1:nu-1)-q2(2:nu))
end subroutine impl
subroutine expl
a=-a
c=-c
b=-a-c
do i=2,nu-1
q2(i)=a(i)*q1(i-1)+b(i)*q1(i)+c(i)*q1(i+1)+r(i)
end do
q2(1)=b(1)*q1(1)+c(1)*q1(2)+r(1)
q2(nu)=b(nu)*q1(nu)+a(nu-1)*q1(nu-1)+r(nu)
fq(2:nu)=kq(2:nu)*rdzu(2:nu)*(q1(1:nu-1)-q1(2:nu))
end subroutine expl
end subroutine diffuse_and_stepx
subroutine diffuse_and_step(q2,q1,fq,dqdt,kq) 13,3
use kvar
; use kcnst
implicit none
real(mrl), dimension(:), intent(inout) :: q2,q1,fq,dqdt,kq
real(mrl), allocatable, target, dimension(:) :: a,b,c,r
allocate(a(nu),b(nu),c(nu),r(nu))
do i=2,nu
a(i)=-kq(i)*dt*rdzu(i)*rdmf(i)
end do
do i=1,nu-1
c(i)=-kq(i+1)*dt*rdzu(i+1)*rdmf(i)
end do
c(nu)=0.
a(1)=0.
b=1-a-c
r=q1*dmfold*rdmf+dt*dqdt
call impl
!call expl
deallocate(a,b,c,r)
contains
subroutine impl 2,2
call tridag
(a,b,c,r,q2)
fq(2:nu)=kq(2:nu)*rdzu(2:nu)*(q2(1:nu-1)-q2(2:nu))
end subroutine impl
subroutine expl
a=-a
c=-c
b=-a-c
do i=2,nu-1
q2(i)=a(i)*q1(i-1)+b(i)*q1(i)+c(i)*q1(i+1)+r(i)
end do
q2(1)=b(1)*q1(1)+c(1)*q1(2)+r(1)
q2(nu)=b(nu)*q1(nu)+a(nu-1)*q1(nu-1)+r(nu)
fq(2:nu)=kq(2:nu)*rdzu(2:nu)*(q1(1:nu-1)-q1(2:nu))
end subroutine expl
end subroutine diffuse_and_step
subroutine prog_tke 1,3
use kvar
; use kcnst
real(mrl), allocatable, target,dimension(:) ::&
a,b,c,r,dudz,dvdz,tmp,elupa,eldna,dummy
real(mrl) fthl_coef,fqt_coef,qvf,qrf,qcf,thf,tmf,exf,ladry,lawet,brmax
integer ibrmax
real(mrl) delzx,elup,sumup,dsumx,eldn,sumdn,thk,thkp,thkm,hitee
allocate(a(nf),b(nf),c(nf),r(nf),dudz(nf),dvdz(nf),tmp(nf),elupa(nu),eldna(nu),dummy(nf))
!hite=get_hite(tke1)
!hite=get_hite2()
!hite=max(2*hite,100.0_mrl)
brmax=0.
dummy=0.
if (time.gt.0.) dummy=buoyflux
do i=2,nf
a(i)=ke(i-1)*dt*rdzf(i-1)*rdmu(i)
end do
do i=1,nf-1
c(i)=ke(i)*dt*rdzf(i)*rdmu(i)
end do
c(nf)=0.
a(1)=0.
do i=2,nu
dudz(i)=rdzu(i)*(u(i)-u(i-1))
dvdz(i)=rdzu(i)*(v(i)-v(i-1))
end do
!do i=1,nu
! debugu2(i)=qv(i)+qc(i)+qr(i)
! debugu3(i)=th(i)-elcp*(qc(i)+qr(i))/exner(i)
!end do
!debugf3=0.
do i=2,nu ! new code
qvf=.5*(qv(i)+qv(i-1))
qcf=.5*(qc(i)+qc(i-1))
qrf=.5*(qr(i)+qr(i-1))
thf=.5*(th(i)+th(i-1))
tmf=.5*(exner(i)*th(i)+exner(i-1)*th(i-1))
exf=.5*(exner(i)+exner(i-1))
! if ((qc(i) > 1.e-12 .or. qc(i-1) > 1.e-12 .or. qr(i) > 1.e12 .or. qr(i-1) > 1.e12 ) .and. wetturb > 0) then !wet
!!!! buoyflux(i)=fth(i)
fthl_coef=(1-(qvf+qrf+qcf)+1.61*qvf*(1+5418./tmf))/&
(1+1.35e7*qvf/(tmf*tmf)) !chi=theta/t
fqt_coef=elcp*fthl_coef/exf-thf !chi=theta/t
! fthl_coef=(1-(qvf+qrf+qcf)+1.61*qvf*(1+5418./tmf))/&
! (1+1.35e7*qvf/(thf*tmf)) !chi=1
! fqt_coef=elcp*fthl_coef-thf !chi=1
! bruntwet(i)=fthl_coef*rdzu(i)*&
! (th(i)-elcp*(qc(i)+qr(i))/exf-th(i-1)+elcp*(qc(i-1)+qr(i-1))/exf)&
! +fqt_coef*rdzu(i)*(qv(i)+qc(i)+qr(i)-qv(i-1)-qc(i-1)-qr(i-1))
bruntwet(i)=fthl_coef*rdzu(i)*&
(thn(i)-elcp*(qcn(i)+qrn(i))/exf-thn(i-1)+elcp*(qcn(i-1)+qrn(i-1))/exf)&
+fqt_coef*rdzu(i)*(qvn(i)+qcn(i)+qrn(i)-qvn(i-1)-qcn(i-1)-qrn(i-1))
buoyfluxwet(i)=-kh(i)*bruntwet(i)
bruntwet(i)=g*bruntwet(i)/th(i)
! buoyflux(i)=-fthl_coef*kh(i)*rdzu(i)*&
! (th(i)-elcp*(qc(i)+qr(i))/exf-th(i-1)+elcp*(qc(i-1)+qr(i-1))/exf)&
! -fqt_coef*kh(i)*rdzu(i)*(qv(i)+qc(i)+qr(i)-qv(i-1)-qc(i-1)-qr(i-1))
! debugf1(i)=-fthl_coef*kh(i)*rdzu(i)*&
! (th(i)-elcp*(qc(i)+qr(i))/exf-th(i-1)+elcp*(qc(i-1)+qr(i-1))/exf)
! debugf2(i)=&
! -fqt_coef*kh(i)*rdzu(i)*(qv(i)+qc(i)+qr(i)-qv(i-1)-qc(i-1)-qr(i-1))
! debugf3(i)=1.
! else !dry
!!!! buoyflux(i)=fth(i)
! bruntdry(i)=rdzu(i)*(1+.61*qvf-qrf-qcf)*(th(i)-th(i-1))&
! +thf*rdzu(i)*(.61*qv(i)-qc(i)-qr(i)-.61*qv(i-1)+qc(i-1)+qr(i-1))
bruntdry(i)=rdzu(i)*(1+.61*qvf-qrf-qcf)*(thn(i)-thn(i-1))&
+thf*rdzu(i)*(.61*qvn(i)-qcn(i)-qrn(i)-.61*qvn(i-1)+qcn(i-1)+qrn(i-1))
if (bruntdry(i)>brmax) then
ibrmax=i
brmax=bruntdry(i)
end if
buoyfluxdry(i)=-kh(i)*bruntdry(i)
if (wetphys.eq.0) buoyfluxdry(i)=fth(i)
if (wetturb.eq.0) buoyfluxdry(i)=fth(i) !new
if (wetturb.eq.0) buoyfluxwet(i)=fth(i) !new
bruntdry(i)=g*bruntdry(i)/th(i)
! buoyflux(i)=fth(i)*(1+.61*qvf-qrf-qcf) &
! +thf*(.61*fqv(i)-fqc(i)-fqr(i))
end do
!Brian 2025:
if (ibrmax.eq.0) ibrmax=1
bruntdry(ibrmax+1)=0.5*(bruntdry(ibrmax)+bruntdry(ibrmax+2))
thetav=th*(1.+.611*qv-qc-qr)
buoyflux(1)=fth(1)*(1+.61*qv(1)-qr(1)-qc(1)) &
+th(1)*(.61*fqv(1))
if (wetturb.eq.0) buoyflux(1)=fth(1) !new
do i=2,nu
brunt(i)=cloudfrac(i)*bruntwet(i)+(1-cloudfrac(i))*bruntdry(i)
buoyflux(i)=cloudfrac(i)*buoyfluxwet(i)+(1-cloudfrac(i))*buoyfluxdry(i)
! debugf1(i)=buoyflux(i)
! debugf2(i)=0.
! debugf3(i)=2.
! end if
end do
! write (69,*) time,buoyflux(10),debugf3(10)
if (mixopt.eq.1) then
!this option sucks
do i=2,nf
! lengthmix(i)=0.25*(1.8*hite*(1-exp(-4*zf(i)/hite)&
! -0.0003*exp(8*zf(i)/hite)))
if (zf(i)/hite.gt.4.) then
lengthmix(i)=1._mrl
else
hitee=hite+40. !desperately trying to get it to entrain
lengthmix(i)=0.25*(1.8*0.9*hitee*(1-exp(-4*zf(i)/(0.9*hitee))&
-0.0003*exp(8*zf(i)/(0.9*hitee))))
end if
! lengthmix(i)=5._mrl
lengthmix(i)=max(lengthmix(i),1._mrl)
lengthdiss(i)=max(lengthmix(i),1._mrl)
!print *,i,lengthdiss(i),lengthmix(i),hite
end do
else if (mixopt.eq.3) then
do k=1,nu-1
! sumup=tke(k)*thetav(k)/g ! wrong mistake found August 29,2001
! sumdn=tke(k)*thetav(k)/g !wrong
sumup=.5*(tke(k)+tke(k+1))*thetav(k)/g
sumdn=.5*(tke(k)+tke(k+1))*thetav(k)/g
elup=0.
eldn=0.
j=k
do while (sumup.gt.0._mrl.and.j.lt.nu)
!thk=thetav(k)
!thkp=thetav(k)
thk =thv(th(k),qv(k),qc(k),qr(k),exner(j))
thkp=thv(th(k),qv(k),qc(k),qr(k),exner(j+1))
delzx=zu(j+1)-zu(j)
dsumx=delzx*.5*(thetav(j)+thetav(j+1)-thkp-thk)
if (sumup-dsumx.gt.0.0.or.thetav(j+1).lt.thetav(j).or.thetav(j).lt.thk-.001) then
sumup=sumup-dsumx
elup=elup+delzx
else
elup=elup+delzx*aproot(0.5*(thetav(j+1)-thetav(j)-thkp+thk),(thetav(j)-thk),-sumup/delzx)
!elup=elup+delzx
sumup=0.
end if
j=j+1
end do
j=k
do while (sumdn.gt.0._mrl.and.j.gt.1)
delzx=zu(j)-zu(j-1)
!thk=thetav(k)
!thkm=thetav(k)
thk =thv(th(k),qv(k),qc(k),qr(k),exner(j))
thkm=thv(th(k),qv(k),qc(k),qr(k),exner(j-1))
dsumx=-delzx*.5*(thetav(j)+thetav(j-1)-thk-thkm)
if (sumdn-dsumx.gt.0.0.or.thetav(j-1).gt.thetav(j).or.thetav(j).gt.thk+.001) then
sumdn=sumdn-dsumx
eldn=eldn+delzx
else
eldn=eldn+delzx*aproot(0.5*(thetav(j)-thetav(j-1)-thk+thkm),-(thetav(j)-thk),-sumdn/delzx)
!eldn=eldn+delzx
sumdn=0.
end if
j=j-1
end do
elupa(k)=elup
eldna(k)=eldn
debugu3(k)=k-j
debugu1(k)=elup
debugu2(k)=eldn
end do
do i=2,nf-1
if (addlength.eq.3) then
lengthmix(i)=.5*sqrt((elupa(i)+elupa(i-1))*(eldna(i)+eldna(i-1)))
else if (addlength.eq.2) then
lengthmix(i)=.5*((elupa(i)+elupa(i-1))+(eldna(i)+eldna(i-1)))
else
lengthmix(i)=.5*((elupa(i)+elupa(i-1))*(eldna(i)+eldna(i-1))) &
/((elupa(i)+elupa(i-1))+(eldna(i)+eldna(i-1)))
end if
lengthdiss(i)=lengthmix(i)
end do
else if (mixopt.eq.2) then
do i=2,nf
lengthdiss(i)=1/(.4*zf(i))+1./300.
lengthmix(i)=lengthdiss(i)
ladry=0.
if (bruntdry(i).gt.0.) ladry=1./(dutchfac*sqrt(tke(i)/bruntdry(i)))
lawet=0.
if (bruntwet(i).gt.0.) lawet=1./(dutchfac*sqrt(tke(i)/bruntwet(i)))
lengthmix(i)=lengthmix(i)+ladry*(1-cloudfrac(i))+lawet*cloudfrac(i)
!if (brunt(i).gt.0.) lengthmix(i)=lengthmix(i)+1./(dutchfac*sqrt(tke(i)/brunt(i)))
lengthdiss(i)=1./lengthdiss(i)
lengthmix(i)=1./lengthmix(i)
lengthmix(i)=max(lengthmix(i),1._mrl)
end do
else
print *, 'mixopt =', mixopt,' not implemented'
stop
end if
lengthmix(1)=25.
lengthdiss(1)=25.
lengthmix(nf)=5.
lengthdiss(nf)=5.
! if (mixopt.eq.2) then
! do k=1,nu
! bruntwet(k)=.5*lengthmix(k)+.5*lengthmix(k+1)
! end do
! do k=2,nf-1
! lengthmix(k)=.5*bruntwet(k)+.5*bruntwet(k-1)
! end do
! end if
tmp(1)=1.
do i=1,nf
!! tmp(i)=(1+dt*sqrt(tke(i))/lengthmix(i))
tmp(i)=1+dt*dissipcoef*sqrt(tke(i))/lengthdiss(i)
end do
buoyflux=1.*buoyflux+.0*dummy
do i=1,nu
r(i)=tke(i)+dt*min(.5*tke(i)/dt, &
max(dtkedt(i)&
-dudz(i)*fu(i)&
-dvdz(i)*fv(i)&
+g*buoyflux(i)/th(i),-.25*tke(i)/dt)& ! new code
)
if (i>1) a(i)=-ke(i-1)*dt*rdzf(i-1)*rdmu(i)
c(i)=-ke(i)*dt*rdzf(i)*rdmu(i)
! a(i)=-.001*dt
! c(i)=-.001*dt
end do
debugf2=dudz
debugf3=fu
a(1)=0.
!c(1)=0.
!c(1)=-ke(1)*dt*rdzf(1)*rdmu(1)
a(nf)=0.
c(nf)=0.
r(1)=tke(1)
r(nf)=tke(nf)
!a=-a
!c=-c
b=-a-c+tmp
call tridag
(a,b,c,r,tken)
do i=1,nf
tken(i)=max(tken(i),0.000001_mrl)
end do
!!TEMP add bullshit filter
!do i=1,nu
!a(i)=.5*(tken(i)+tken(i+1))
!end do
!do i=2,nf-1
!tken(i)=.5*(a(i-1)+a(i))
!end do
!print *, ' tke=',tken(2)
deallocate(a,b,c,r,dudz,dvdz,tmp,elupa,eldna,dummy)
end subroutine prog_tke
function thv(th,qvv,qc,qr,exner),4
use kvar
, only: mrl,psur,n,wetphys, elcp; use kcnst
; use micromod
real(mrl) th,ql,exner,dqv,t,p,esl,desdtl,qsat,qc,qr,qv,thv,qvv
qv=qvv
if (wetphys.eq.0.or.1.eq.0) then
thv=th*(1+.611*qv-qr-qc)
else
ql=qc+qr
t=th*exner
p=psur*(exner**3.5)
do n=1,1
call sat_water
(t,p,esl,desdtl,qsat)
dqv=(qsat-qv)/(1+.611*desdtl*elcp/p)
dqv=min(ql,dqv)
qv=qv+dqv
ql=ql-dqv
t=t-elcp*dqv
end do
thv=(t/exner)*(1.+.611*qv-ql)
end if
return
end function thv
function aproot(a,b,c),2
use kvar
; use kcnst
real(mrl) a,b,c,aproot
if (abs(a).lt.0.001*abs(b)) then
!print *, a,b,c,aproot
aproot=-c/b
else
aproot=(-b+sqrt(b*b-4*a*c))/(2*a)
end if
end function aproot
function f_get_hite(dum) ,2
!implicit none
use kvar
; use kcnst
real(mrl), dimension(:), intent(in) :: dum
real(mrl) f_get_hite
m=1
do while (dum(m).gt.0.5*dum(1))
m=m+1
end do
f_get_hite=zu(m-1)+(0.5*dum(1)-dum(m-1))/(dum(m)-dum(m-1))*(zu(m)-zu(m-1))
return
end function f_get_hite
function get_hite2() ,2
!implicit none
use kvar
; use kcnst
real(mrl) grad,dummy,get_hite2
grad=0.
do n=1,nu-1
!dummy=(th1(n+1)-th1(n))/(zu(n+1)-zu(n))
dummy=(smoke(n+1)-smoke(n))/(zu(n+1)-zu(n))
!if (dummy.gt.grad) then
if (dummy.lt.grad) then
get_hite2=zf(n+1)
grad=dummy
end if
end do
return
end function get_hite2
function f_get_hite3(b,ifilt) 1,3
use kvar
; use kcnst
implicit none
real (mrl), intent(in) :: b(:)
real (mrl), allocatable :: dummy(:),grad(:),fummy(:)
real(mrl) extrem,f_get_hite3
integer iex,ifilt,nn
nn=size(b)
allocate(dummy(nn),fummy(nn),grad(nn+1))
dummy=b
do k=1,ifilt
fummy(1)=.5*(dummy(1)+dummy(2))
fummy(nn)=.5*(dummy(nn)+dummy(nn-1))
do n=2,nn-1
fummy(n)=.5*dummy(n)+.25*dummy(n-1)+.25*dummy(n+1)
end do
dummy=fummy
end do
grad(1)=0.
grad(nn+1)=0.
do n=2,nn
grad(n)=rdzu(n)*(dummy(n)-dummy(n-1))
end do
extrem=grad(1)
iex=1
do n=2,nu
if (grad(n).lt.extrem) then
iex=n
extrem=grad(n)
end if
end do
!Brian 2025:
if (iex.lt.2) iex=2
f_get_hite3=f_paramin
(zf(iex-1),zf(iex),zf(iex+1),&
grad(iex-1),grad(iex),grad(iex+1))
deallocate(dummy,grad)
return
end function f_get_hite3
subroutine tke_to_rho_k 1,2
use kvar
; use kcnst
if (khval.le.0.) then
do i=2,nu
kh(i)=diffucoef*sqrt(2*tke(i))*lengthmix(i)*(rho(i-1)+rho(i))*.5
end do
else
kh=khval
end if
kh(1)=0.
kh(nf)=0.
km=prandtl*kh
do i=1,nu
ke(i)=.5*(km(i)+km(i+1))
end do
return
end subroutine tke_to_rho_k
subroutine repoint 4,2
use kvar
; use kcnst
ipnt=mod(istep,2)+1
jpnt=mod(istep+1,2)+1
include '../inc/repoint.inc'
return
end subroutine repoint
subroutine tridag(a,b,c,r,u) 3,1
use precis
implicit none
real(mrl), dimension(:), intent(in) :: a,b,c,r
real(mrl), dimension(:), intent(out) :: u
real(mrl), dimension(size(b)) :: gam
integer :: n,j
real(mrl) :: bet
n=size(b)
bet=b(1)
if (bet == 0.0) print *,'error at code stage 1'
u(1)=r(1)/bet
do j=2,n
gam(j)=c(j-1)/bet
!bet=b(j)-a(j-1)*gam(j) ! bloody bastards
bet=b(j)-a(j)*gam(j)
if (bet == 0.0) &
print *,' error at code stage 2'
!u(j)=(r(j)-a(j-1)*u(j-1))/bet !bloody bastards
u(j)=(r(j)-a(j)*u(j-1))/bet
end do
do j=n-1,1,-1
u(j)=u(j)-gam(j+1)*u(j+1)
end do
end subroutine tridag
function f_paramin(a,b,c,fa,fb,fc) 1,1
use precis
real(mrl) a,b,c,fa,fb,fc,f_paramin
f_paramin=b-0.5*((b-a)**2*(fb-fc)-(b-c)**2*(fb-fa))/&
((b-a)*(fb-fc)-(b-c)*(fb-fa))
return
end function f_paramin
subroutine update_time 1,3
use kvar
; use kcnst
time=time+dt
istep=istep+1
call repoint
include '../inc/surface.inc'
end subroutine update_time
end module kphysmod