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

move setting of flag from run to init phase #2

Merged
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
26 changes: 14 additions & 12 deletions physics/cu_c3_driver.F90
Original file line number Diff line number Diff line change
Expand Up @@ -30,21 +30,31 @@ module cu_c3_driver
!! \htmlinclude cu_c3_driver_init.html
!!
subroutine cu_c3_driver_init(imfshalcnv, imfshalcnv_c3, imfdeepcnv, &
imfdeepcnv_c3,mpirank, mpiroot, errmsg, errflg)
imfdeepcnv_c3,progsigma, cnx, mpirank, mpiroot, &
errmsg, errflg)

implicit none

integer, intent(in) :: imfshalcnv, imfshalcnv_c3
integer, intent(in) :: imfdeepcnv, imfdeepcnv_c3
integer, intent(in) :: mpirank
integer, intent(in) :: mpiroot
integer, intent(in) :: cnx
logical, intent(inout) :: progsigma
character(len=*), intent( out) :: errmsg
integer, intent( out) :: errflg

! initialize ccpp error handling variables
errmsg = ''
errflg = 0

if(progsigma)then
if(cnx < 384)then
progsigma=.false.
write(*,*)'Forcing prognostic closure to .false. due to coarse resolution'
endif
endif

end subroutine cu_c3_driver_init

!
Expand All @@ -58,7 +68,7 @@ end subroutine cu_c3_driver_init
!!
!>\section gen_c3_driver Grell-Freitas Cumulus Scheme Driver General Algorithm
subroutine cu_c3_driver_run(ntracer,garea,im,km,dt,flag_init,flag_restart,&
do_ca,progsigma,cnx,cactiv,cactiv_m,g,cp,fv,r_d,xlv,r_v,forcet, &
do_ca,progsigma,cactiv,cactiv_m,g,cp,fv,r_d,xlv,r_v,forcet, &
forceqv_spechum,phil,delp,raincv,tmf,qmicro,sigmain, &
betascu,betamcu,betadcu,qv_spechum,t,cld1d,us,vs,t2di,w, &
qv2di_spechum,p2di,psuri, &
Expand Down Expand Up @@ -93,14 +103,14 @@ subroutine cu_c3_driver_run(ntracer,garea,im,km,dt,flag_init,flag_restart,&
integer :: ishallow_g3 ! depend on imfshalcnv
!-------------------------------------------------------------
integer :: its,ite, jts,jte, kts,kte
integer, intent(in ) :: im,km,ntracer,cnx
integer, intent(in ) :: im,km,ntracer
integer, intent(in ) :: ichoice_in,ichoicem_in,ichoice_s_in
logical, intent(in ) :: flag_init, flag_restart, do_mynnedmf
logical, intent(in ) :: flag_for_scnv_generic_tend,flag_for_dcnv_generic_tend, &
do_ca
real (kind=kind_phys), intent(in) :: g,cp,fv,r_d,xlv,r_v,betascu,betamcu,betadcu
logical, intent(in ) :: ldiag3d
logical, intent(inout) :: progsigma
logical, intent(in ) :: progsigma
real(kind=kind_phys), intent(inout) :: dtend(:,:,:)
!$acc declare copy(dtend)
integer, intent(in) :: dtidx(:,:), &
Expand Down Expand Up @@ -280,14 +290,6 @@ subroutine cu_c3_driver_run(ntracer,garea,im,km,dt,flag_init,flag_restart,&
!$acc end kernels
endif


if(progsigma)then
if(cnx < 384)then
progsigma=.false.
write(*,*)'Forcing prognostic closure to .false. due to coarse resolution'
endif
endif

if(ldiag3d) then
if(flag_for_dcnv_generic_tend) then
cliw_deep_idx=0
Expand Down
21 changes: 14 additions & 7 deletions physics/cu_c3_driver.meta
Original file line number Diff line number Diff line change
Expand Up @@ -49,6 +49,20 @@
dimensions = ()
type = integer
intent = in
[progsigma]
standard_name = do_prognostic_updraft_area_fraction
long_name = flag for prognostic sigma in cumuls scheme
units = flag
dimensions = ()
type = logical
intent = inout
[cnx]
standard_name = number_of_x_points_for_current_cubed_sphere_tile
long_name = number of points in x direction for this cubed sphere face
units = count
dimensions = ()
type = integer
intent = in
[errmsg]
standard_name = ccpp_error_message
long_name = error message for error handling in CCPP
Expand Down Expand Up @@ -133,13 +147,6 @@
units = flag
dimensions = ()
type = logical
intent = inout
[cnx]
standard_name = number_of_x_points_for_current_cubed_sphere_tile
long_name = number of points in x direction for this cubed sphere face
units = count
dimensions = ()
type = integer
intent = in
[cactiv]
standard_name = counter_for_grell_freitas_convection
Expand Down
Loading