Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

D80R scheme #63

Open
wants to merge 1 commit into
base: to4.4_Fredrik
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
34 changes: 21 additions & 13 deletions src/modbudget.f90
Original file line number Diff line number Diff line change
Expand Up @@ -33,7 +33,7 @@ module modbudget
PUBLIC :: initbudget,budgetstat,exitbudget
save
!NetCDF variables
integer,parameter :: nvar = 18
integer,parameter :: nvar = 19
character(80),dimension(nvar,4) :: ncname

real :: dtav, timeav
Expand Down Expand Up @@ -67,6 +67,7 @@ module modbudget
real, allocatable :: sbtkeav(:) ! Must be module global to calculate dE/dt in combination with tkeb
real, allocatable :: ekmmn(:) !< Turbulent exchange coefficient momentum
real, allocatable :: khkmmn(:) !< Kh / Km, in post-processing used to determine filter-grid ratio
real, allocatable :: zltmn(:) !< length scale parameter


logical :: ltkeb !Switch to tell if the tke at beg of av periode has been stored
Expand Down Expand Up @@ -128,7 +129,7 @@ subroutine initbudget
!time averaged fields, subgrid TKE
allocate(sbtkemn(k1),sbtkeb(k1),sbtkeav(k1),sbshrmn(k1),sbbuomn(k1))
allocate(sbstormn(k1),sbbudgmn(k1),sbresidmn(k1),&
sbdissmn(k1),ekmmn(k1),khkmmn(k1))
sbdissmn(k1),ekmmn(k1),khkmmn(k1),zltmn(k1))


!Setting time mean variables to zero
Expand All @@ -137,7 +138,7 @@ subroutine initbudget
tkeb=0. !tke at start of averaging period
sbtkemn=0.;sbshrmn=0.;sbbuomn=0.;sbstormn=0.;
sbbudgmn=0.;sbresidmn=0.;sbdissmn=0.;
ekmmn=0.;khkmmn=0.
ekmmn=0.;khkmmn=0.;zltmn=0.;
sbtkeb=0.;
ltkeb=.false. ; lsbtkeb=.false.

Expand Down Expand Up @@ -173,6 +174,7 @@ subroutine initbudget
call ncinfo(ncname(16,:),'sbresid','Subgrid Residual = budget - storage','m/s^2','tt')
call ncinfo(ncname(17,:),'ekm','Turbulent exchange coefficient momentum','m/s^2','tt')
call ncinfo(ncname(18,:),'khkm ','Kh / Km, in post-processing used to determine filter-grid ratio','m/s^2','tt')
call ncinfo(ncname(19,:),'zlt ','Length scale parameter','m','tt')
call define_nc( ncid_prof, NVar, ncname)
end if

Expand Down Expand Up @@ -668,27 +670,27 @@ end subroutine do_genbudget
!> Performs the SFS - budget calculations
subroutine do_gensbbudget
use modglobal, only : i1,j1,ih,jh,k1,ijtot
use modsubgriddata, only : ekm,ekh,sbdiss,sbshr,sbbuo
use modsubgriddata, only : ekm,ekh,zlt,sbdiss,sbshr,sbbuo
use modfields, only : e120,rhobf
use modmpi, only : slabsum,comm3d,my_real, mpi_sum,mpierr
!----------------------------
! 1.1 Declare allocatable
!----------------------------
integer i,j,k
real, allocatable :: sbshrav(:),sbbuoav(:) ,sbdissav(:),&
ekmav(:),khkmav(:)
ekmav(:),khkmav(:),zltav(:)
real, allocatable :: sbtkeavl(:),khkmavl(:)

!----------------------------
! 1.2 Allocate variables
!----------------------------
allocate(sbshrav(k1),sbbuoav(k1),sbdissav(k1),ekmav(k1),khkmav(k1))
allocate(sbshrav(k1),sbbuoav(k1),sbdissav(k1),ekmav(k1),khkmav(k1),zltav(k1))
allocate(sbtkeavl(k1),khkmavl(k1))

!----------------------------
! 2. Reset variables
!----------------------------
sbtkeav=0.;sbshrav=0.;sbbuoav=0.;sbdissav=0;ekmav=0.;khkmav=0.
sbtkeav=0.;sbshrav=0.;sbbuoav=0.;sbdissav=0;ekmav=0.;khkmav=0.;zltav=0.
sbtkeavl=0.;khkmavl=0.

!----------------------------
Expand All @@ -710,13 +712,15 @@ subroutine do_gensbbudget
call slabsum(sbbuoav ,1,k1,sbbuo ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(sbdissav,1,k1,sbdiss,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(ekmav ,1,k1,ekm ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(zltav ,1,k1,zlt ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call MPI_ALLREDUCE(khkmavl, khkmav, k1, MY_REAL, MPI_SUM, comm3d,mpierr)
call MPI_ALLREDUCE(sbtkeavl,sbtkeav,k1, MY_REAL, MPI_SUM, comm3d,mpierr)
sbshrav = sbshrav / ijtot
sbbuoav = sbbuoav / ijtot
sbdissav = sbdissav/ ijtot
sbtkeav = sbtkeav / ijtot
ekmav = ekmav / ijtot
zltav = zltav / ijtot
khkmav = khkmav / ijtot
!storage: tke at beginning of averaging period
if (.not.(lsbtkeb)) then
Expand All @@ -735,10 +739,11 @@ subroutine do_gensbbudget
sbtkemn(k) = sbtkemn(k) + rhobf(k)*sbtkeav(k)
ekmmn(k) = ekmmn(k) + ekmav(k)
khkmmn(k) = khkmmn(k) + khkmav(k)
zltmn(k) = zltmn(k) + zltav(k)
enddo

!Deallocate
deallocate(sbshrav,sbbuoav,sbdissav,ekmav,khkmav)
deallocate(sbshrav,sbbuoav,sbdissav,ekmav,khkmav,zltav)
deallocate(sbtkeavl,khkmavl)
end subroutine do_gensbbudget

Expand Down Expand Up @@ -773,6 +778,7 @@ subroutine writebudget
sbbuomn = sbbuomn /nsamples
sbdissmn = sbdissmn /nsamples
ekmmn = ekmmn /nsamples
zltmn = zltmn /nsamples
khkmmn = khkmmn /nsamples
sbstormn = (sbtkeav-sbtkeb)/timeav
sbbudgmn = sbshrmn+sbbuomn+sbdissmn
Expand Down Expand Up @@ -821,13 +827,13 @@ subroutine writebudget
write (ifoutput,'(A/2A/3A)') &
'#---------------------------------------------------------------' &
,'#LEV HEIGHT | SBTKE SBSHEAR BUOYANCY SBDISS' &
,' SBSTORAGE SBBUDGET SBRESID EKM KH/KM '&
,' SBSTORAGE SBBUDGET SBRESID EKM KH/KM LAMBDA'&
,'# (M) | ' &
,'(-------------------------------------------------- (M/S)^2 -------------' &
,'-----------------------------------)'
,'-----------------------------------------)'


write(ifoutput,'(I3,F9.3,9E12.4)') &
write(ifoutput,'(I3,F9.3,10E12.4)') &
(k, &
zf (k), &
sbtkemn (k), &
Expand All @@ -839,6 +845,7 @@ subroutine writebudget
sbresidmn (k), &
ekmmn (k), & !!!
khkmmn (k), &
zltmn (k), &
k=1,kmax)
close(ifoutput)
endif !endif myid==0
Expand All @@ -861,6 +868,7 @@ subroutine writebudget
vars(:,16) =sbresidmn
vars(:,17) =ekmmn
vars(:,18) =khkmmn
vars(:,19) =zltmn
call writestat_nc(ncid_prof,nvar,ncname,vars(1:kmax,:),nrec_prof,kmax)
end if
!Reset time mean variables; resolved TKE
Expand All @@ -869,7 +877,7 @@ subroutine writebudget
ltkeb=.false.
!Reset time mean variables; subgrid TKE
sbtkemn=0.;sbtkeb=0.;sbshrmn=0.;sbbuomn=0.;sbdissmn=0.;sbtkeb=0.
sbstormn=0.;sbbudgmn=0.;sbresidmn=0.;ekmmn=0.;khkmmn=0.;
sbstormn=0.;sbbudgmn=0.;sbresidmn=0.;ekmmn=0.;khkmmn=0.;zltmn=0.;
lsbtkeb=.false.
end subroutine writebudget

Expand All @@ -882,7 +890,7 @@ subroutine exitbudget

deallocate(tkemn,tkeb,tkeav,shrmn,buomn,trspmn,ptrspmn,stormn,budgmn,residmn,dissmn)
deallocate(sbtkemn,sbshrmn,sbbuomn,sbstormn,sbbudgmn,sbresidmn,sbdissmn,sbtkeb,sbtkeav)
deallocate(ekmmn,khkmmn)
deallocate(ekmmn,khkmmn,zltmn)
end subroutine exitbudget

end module modbudget
Expand Down
9 changes: 6 additions & 3 deletions src/modsubgrid.f90
Original file line number Diff line number Diff line change
Expand Up @@ -36,7 +36,7 @@ module modsubgrid

contains
subroutine initsubgrid
use modglobal, only : ih,i1,jh,j1,k1,delta,zf,fkar,pi
use modglobal, only : ih,i1,jh,j1,k1,delta,zh,zf,fkar,pi
use modmpi, only : myid

implicit none
Expand Down Expand Up @@ -106,7 +106,7 @@ subroutine subgridnamelist
integer :: ierr

namelist/NAMSUBGRID/ &
ldelta,lmason, cf,cn,Rigc,Prandtl,lsmagorinsky,cs,nmason,ch1,ch2,cm,ce1,ce2
ldelta,lmason,lD80R, cf,cn,Rigc,Prandtl,lsmagorinsky,cs,nmason,ch1,ch2,cm,ce1,ce2

if(myid==0)then
open(ifnamopt,file=fname_options,status='old',iostat=ierr)
Expand All @@ -122,6 +122,7 @@ subroutine subgridnamelist

call MPI_BCAST(ldelta ,1,MPI_LOGICAL,0,comm3d,mpierr)
call MPI_BCAST(lmason ,1,MPI_LOGICAL,0,comm3d,mpierr)
call MPI_BCAST(lD80R ,1,MPI_LOGICAL,0,comm3d,mpierr)
call MPI_BCAST(nmason ,1,MY_REAL ,0,comm3d,mpierr)
call MPI_BCAST(lsmagorinsky,1,MPI_LOGICAL,0,comm3d,mpierr)
call MPI_BCAST(cs ,1,MY_REAL ,0,comm3d,mpierr)
Expand Down Expand Up @@ -200,7 +201,7 @@ subroutine closure
! |
!-----------------------------------------------------------------|

use modglobal, only : i1,j1,kmax,k1,ih,jh,i2,j2,delta,ekmin,grav,zf,fkar,deltai, &
use modglobal, only : i1,j1,kmax,k1,ih,jh,i2,j2,delta,ekmin,grav,zh,zf,fkar,deltai, &
dxi,dyi,dzf,dzh
use modfields, only : dthvdz,e120,u0,v0,w0,thvf
use modsurfdata, only : dudz,dvdz,z0m
Expand Down Expand Up @@ -316,6 +317,8 @@ subroutine closure

if (lmason) zlt(i,j,k) = (1. / zlt(i,j,k) ** nmason + 1. / ( fkar * (zf(k) + z0m(i,j)))**nmason) ** (-1./nmason)

if (lD80R) zlt(i,j,k) = 1/(1/(0.4*zh(k))+1/(cn*e120(i,j,k)/sqrt(grav/thvf(k)*abs(dthvdz(i,j,k)))))

ekm(i,j,k) = cm * zlt(i,j,k) * e120(i,j,k)
ekh(i,j,k) = (ch1 + ch2 * zlt(i,j,k)*deltai(k)) * ekm(i,j,k)

Expand Down
1 change: 1 addition & 0 deletions src/modsubgriddata.f90
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,7 @@ module modsubgriddata
! private
logical :: ldelta = .false. !< switch for subgrid length formulation (on/off)
logical :: lmason = .false. !< switch for decreased length scale near the surface
logical :: lD80R = .false. !< switch for D80-R scheme
logical :: lsmagorinsky = .false. !< switch for smagorinsky subgrid scheme
real :: cf = 2.5 !< filter constant
real :: Rigc = 0.25 !< critical Richardson number
Expand Down
2 changes: 1 addition & 1 deletion src/modtimedep.f90
Original file line number Diff line number Diff line change
Expand Up @@ -432,7 +432,7 @@ subroutine timedepsurf
do while(rtimee>timeflux(t))
t=t+1
end do
if (rtimee>timeflux(t)) then
if (rtimee>timeflux(t-1)) then
t=t-1
endif

Expand Down
1 change: 1 addition & 0 deletions tags
Original file line number Diff line number Diff line change
Expand Up @@ -1740,6 +1740,7 @@ llimit src/modradfull.f90 /^ real :: llimit,/;" k type:band_
llimit src/modradfull.f90 /^ elemental real function llimit(/;" f module:modradfull
llsadv src/modglobal.f90 /^ logical :: llsadv /;" v module:modglobal
lmason src/modsubgriddata.f90 /^ logical :: lmason /;" v module:modsubgriddata
lD80R src/modsubgriddata.f90 /^ logical :: lD80R /;" v module:modsubgriddata
lmicrostat src/modbulkmicrostat.f90 /^ logical :: lmicrostat /;" v module:modbulkmicrostat
lmicrostat src/modsimpleicestat.f90 /^ logical :: lmicrostat /;" v module:modsimpleicestat
lmoist src/modglobal.f90 /^ logical :: lmoist /;" v module:modglobal
Expand Down