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

Handling of scalar variables #136

Merged
merged 2 commits into from
Nov 7, 2024
Merged
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
43 changes: 15 additions & 28 deletions src/modcrosssection.f90
Original file line number Diff line number Diff line change
Expand Up @@ -284,7 +284,7 @@ subroutine wrtvert
use modstat_nc, only : lnetcdf, writestat_nc
implicit none

integer i,k,n
integer i,k,n,isv
character(20) :: name

real, allocatable :: thv0(:,:),vars(:,:,:),buoy(:,:)
Expand Down Expand Up @@ -362,14 +362,10 @@ subroutine wrtvert
vars(:,:,6) = qtm(2:i1,crossplane(cross),1:kmax)
vars(:,:,7) = ql0(2:i1,crossplane(cross),1:kmax)
vars(:,:,8) = buoy(2:i1,1:kmax)
if(nsv>1) then
vars(:,:,9) = svm(2:i1,crossplane(cross),1:kmax,2)
vars(:,:,10) = svm(2:i1,crossplane(cross),1:kmax,1)
else
vars(:,:,9) = 0.
vars(:,:,10) = 0.
end if
vars(:,:,11) = e120(2:i1,crossplane(cross),1:kmax)
vars(:,:,9) = e120(2:i1,crossplane(cross),1:kmax)
do isv = 1,nsv
vars(:,:,9+isv) = svm(2:i1,crossplane(cross),1:kmax,isv)
end do
call writestat_nc(ncid1(cross),1,tncname1,(/rtimee/),nrec1(cross),.true.)
call writestat_nc(ncid1(cross),nvar,ncname1(1:nvar,:),vars,nrec1(cross),imax,kmax)
end if
Expand All @@ -386,12 +382,11 @@ subroutine wrthorz
use modfields, only : um,vm,wm,thlm,qtm,svm,thl0,qt0,ql0,e120,exnf,thvf
use modmpi, only : cmyid
use modstat_nc, only : lnetcdf, writestat_nc
use modmicrodata, only : iqr,inr
implicit none


! LOCAL
integer i,j,n
integer i,j,n,isv
character(40) :: name
real, allocatable :: thv0(:,:,:),vars(:,:,:),buoy(:,:,:)

Expand Down Expand Up @@ -468,14 +463,10 @@ subroutine wrthorz
vars(:,:,6) = qtm(2:i1,2:j1,crossheight(cross))
vars(:,:,7) = ql0(2:i1,2:j1,crossheight(cross))
vars(:,:,8) = buoy(2:i1,2:j1,cross)
if(nsv>1) then
vars(:,:,9) = svm(2:i1,2:j1,crossheight(cross),iqr)
vars(:,:,10) = svm(2:i1,2:j1,crossheight(cross),inr)
else
vars(:,:,9) = 0.
vars(:,:,10) = 0.
end if
vars(:,:,11) = e120(2:i1,2:j1,crossheight(cross))
vars(:,:,9) = e120(2:i1,2:j1,crossheight(cross))
do isv = 1,nsv
vars(:,:,9+isv) = svm(2:i1,2:j1,crossheight(cross),isv)
end do
call writestat_nc(ncid2(cross),1,tncname2,(/rtimee/),nrec2(cross),.true.)
call writestat_nc(ncid2(cross),nvar,ncname2(1:nvar,:),vars,nrec2(cross),imax,jmax)
end do
Expand All @@ -496,7 +487,7 @@ subroutine wrtorth


! LOCAL
integer j,k,n
integer j,k,n,isv
character(21) :: name

real, allocatable :: thv0(:,:),vars(:,:,:),buoy(:,:)
Expand Down Expand Up @@ -581,14 +572,10 @@ subroutine wrtorth
vars(:,:,6) = qtm(crossortho(cross),2:j1,1:kmax)
vars(:,:,7) = ql0(crossortho(cross),2:j1,1:kmax)
vars(:,:,8) = buoy(2:j1,1:kmax)
if(nsv>1) then
vars(:,:,9) = svm(crossortho(cross),2:j1,1:kmax,2)
vars(:,:,10) = svm(crossortho(cross),2:j1,1:kmax,1)
else
vars(:,:,9) = 0.
vars(:,:,10) = 0.
end if
vars(:,:,11) = e120(crossortho(cross),2:j1,1:kmax)
vars(:,:,9) = e120(crossortho(cross),2:j1,1:kmax)
do isv = 1,nsv
vars(:,:,9+isv) = svm(crossortho(cross),2:j1,1:kmax,isv)
end do
call writestat_nc(ncid3(cross),1,tncname3,(/rtimee/),nrec3(cross),.true.)
call writestat_nc(ncid3(cross),nvar,ncname3(1:nvar,:),vars,nrec3(cross),jmax,kmax)
end if
Expand Down
4 changes: 2 additions & 2 deletions src/modgenstat.f90
Original file line number Diff line number Diff line change
Expand Up @@ -744,7 +744,7 @@ subroutine do_genstat
clwav(1) = clwav(1) + ql0(i,j,1) * ilratio
cliav(1) = cliav(1) + ql0(i,j,1) * (1-ilratio)

if (nsv > 1) then
if (iqr > 0) then
ilratio = max(0._field_r,min(1._field_r,(tmp0(i,j,1)-tdnrsg)/(tuprsg-tdnrsg)))
plwav(1) = plwav(1) + sv0(i,j,1,iqr) * ilratio
pliav(1) = pliav(1) + sv0(i,j,1,iqr) * (1-ilratio)
Expand Down Expand Up @@ -850,7 +850,7 @@ subroutine do_genstat
clwav_s = clwav_s + ql0(i,j,k) * ilratio
cliav_s = cliav_s + ql0(i,j,k) * (1-ilratio)

if (nsv > 1) then
if (iqr > 0) then
ilratio = max(0._field_r,min(1._field_r,(tmp0(i,j,k)-tdnrsg)/(tuprsg-tdnrsg)))
plwav_s = plwav_s + sv0(i,j,k,iqr) * ilratio
pliav_s = pliav_s + sv0(i,j,k,iqr) * (1-ilratio)
Expand Down
2 changes: 1 addition & 1 deletion src/modmicrodata.f90
Original file line number Diff line number Diff line change
Expand Up @@ -50,7 +50,7 @@ module modmicrodata

logical :: l_lognormal = .false. !< log param of rain terminal velocities for rain sedim

integer :: inr, iqr
integer :: inr = -1, iqr = -1

real(field_r), parameter :: D0_kk = 50e-6 & !< diameter sep. cloud and prec. in KK00 scheme
,qcmin = 1.0e-7 & !< Cloud specific mixing ratio treshold for calculations
Expand Down
8 changes: 2 additions & 6 deletions src/modmicrophysics.f90
Original file line number Diff line number Diff line change
Expand Up @@ -86,19 +86,15 @@ subroutine initmicrophysics
case(imicro_none)
case(imicro_drizzle)
case(imicro_bulk)
!if (nsv < 2) STOP "ERROR: Bulk microphysics requires nsv >=2"
call initbulkmicro
call initbulkmicro
case(imicro_bin)
! call initbinmicro
case(imicro_sice)
if (nsv < 2) STOP "ERROR: Simple ice microphysics requires nsv >=2"
call initsimpleice
case(imicro_sice2)
if (nsv < 2) STOP "ERROR: Simple ice microphysics requires nsv >=2"
call initsimpleice2
case(imicro_bulk3) !#sb3
if (nsv < 12) STOP "ERROR: Full Seifer-Beheng microphysics requires nsv >=12" !#sb3
call initbulkmicro3 !#sb3
call initbulkmicro3 !#sb3
case(imicro_user)
end select
end subroutine initmicrophysics
Expand Down
Loading