Skip to content
This repository has been archived by the owner on Oct 23, 2020. It is now read-only.

Remove the call to the cloud microphysics driver from the dynamics RK routine #1410

Open
wants to merge 2 commits into
base: atmosphere/develop
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
68 changes: 0 additions & 68 deletions src/core_atmosphere/dynamics/mpas_atm_time_integration.F
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,6 @@ module atm_time_integration
use mpas_timer

#ifdef DO_PHYSICS
use mpas_atmphys_driver_microphysics
use mpas_atmphys_todynamics
use mpas_atmphys_utilities
#endif
Expand Down Expand Up @@ -165,8 +164,6 @@ subroutine atm_srk3(domain, dt, itimestep)
logical, pointer :: config_positive_definite
logical, pointer :: config_monotonic
real (kind=RKIND), pointer :: config_dt
character (len=StrKIND), pointer :: config_microp_scheme
character (len=StrKIND), pointer :: config_convection_scheme

integer, pointer :: num_scalars, index_qv, nCells, nCellsSolve, nEdges, nEdgesSolve, nVertices, nVerticesSolve, nVertLevels

Expand Down Expand Up @@ -198,8 +195,6 @@ subroutine atm_srk3(domain, dt, itimestep)
real (kind=RKIND), dimension(:,:), pointer :: u, uReconstructZonal, uReconstructMeridional, uReconstructX, uReconstructY, uReconstructZ
real (kind=RKIND), dimension(:,:,:), pointer :: scalars, scalars_1, scalars_2

real (kind=RKIND), dimension(:,:), pointer :: rqvdynten

logical, parameter :: debug = .false.


Expand All @@ -212,8 +207,6 @@ subroutine atm_srk3(domain, dt, itimestep)
call mpas_pool_get_config(domain % blocklist % configs, 'config_positive_definite', config_positive_definite)
call mpas_pool_get_config(domain % blocklist % configs, 'config_monotonic', config_monotonic)
call mpas_pool_get_config(domain % blocklist % configs, 'config_dt', config_dt)
call mpas_pool_get_config(domain % blocklist % configs, 'config_microp_scheme', config_microp_scheme)
call mpas_pool_get_config(domain % blocklist % configs, 'config_convection_scheme', config_convection_scheme)
call mpas_pool_get_config(domain % blocklist % configs, 'config_IAU_option', config_IAU_option)

! config variables for dynamics-transport splitting, WCS 18 November 2014
Expand Down Expand Up @@ -1249,74 +1242,13 @@ subroutine atm_srk3(domain, dt, itimestep)
block => block % next
end do

!
! call to parameterizations of cloud microphysics. calculation of the tendency of water vapor to horizontal and
! vertical advection needed for the Tiedtke parameterization of convection.
!

#ifdef DO_PHYSICS
block => domain % blocklist
do while(associated(block))

call mpas_pool_get_subpool(block % structs, 'mesh', mesh)
call mpas_pool_get_subpool(block % structs, 'state', state)
call mpas_pool_get_subpool(block % structs, 'diag', diag)
call mpas_pool_get_subpool(block % structs, 'diag_physics', diag_physics)
call mpas_pool_get_subpool(block % structs, 'tend_physics', tend_physics)
call mpas_pool_get_subpool(block % structs, 'tend', tend)
call mpas_pool_get_array(state, 'scalars', scalars_1, 1)
call mpas_pool_get_array(state, 'scalars', scalars_2, 2)
call mpas_pool_get_dimension(state, 'index_qv', index_qv)

call mpas_pool_get_dimension(block % dimensions, 'nThreads', nThreads)

call mpas_pool_get_dimension(block % dimensions, 'cellSolveThreadStart', cellSolveThreadStart)
call mpas_pool_get_dimension(block % dimensions, 'cellSolveThreadEnd', cellSolveThreadEnd)

if(config_convection_scheme == 'cu_grell_freitas' .or. &
config_convection_scheme == 'cu_tiedtke' .or. &
config_convection_scheme == 'cu_ntiedtke') then

call mpas_pool_get_array(tend_physics, 'rqvdynten', rqvdynten)

!NOTE: The calculation of the tendency due to horizontal and vertical advection for the water vapor mixing ratio
!requires that the subroutine atm_advance_scalars_mono was called on the third Runge Kutta step, so that a halo
!update for the scalars at time_levs(1) is applied. A halo update for the scalars at time_levs(2) is done above.
if (config_monotonic) then
rqvdynten(:,:) = ( scalars_2(index_qv,:,:) - scalars_1(index_qv,:,:) ) / config_dt
else
rqvdynten(:,:) = 0._RKIND
end if
end if

!simply set to zero negative mixing ratios of different water species (for now):
where ( scalars_2(:,:,:) < 0.0) &
scalars_2(:,:,:) = 0.0

!call microphysics schemes:
if (trim(config_microp_scheme) /= 'off') then
call mpas_timer_start('microphysics')
!$OMP PARALLEL DO
do thread=1,nThreads
call driver_microphysics ( block % configs, mesh, state, 2, diag, diag_physics, tend, itimestep, &
cellSolveThreadStart(thread), cellSolveThreadEnd(thread))
end do
!$OMP END PARALLEL DO
call mpas_timer_stop('microphysics')
end if
block => block % next
end do

!
! Note: A halo exchange for 'exner' here as well as at the end of
! the first (n-1) dynamics subcycles can substitute for the exchange at
! the beginning of each dynamics subcycle. Placing halo exchanges here
! and at the end of dynamics subcycles may in future allow for aggregation
! of the 'exner' exchange with other exchanges.
!
#endif

call summarize_timestep(domain)

end subroutine atm_srk3

Expand Down
28 changes: 25 additions & 3 deletions src/core_atmosphere/mpas_atm_core.F
Original file line number Diff line number Diff line change
Expand Up @@ -825,7 +825,7 @@ subroutine atm_do_timestep(domain, dt, itimestep)
use atm_time_integration
#ifdef DO_PHYSICS
use mpas_atmphys_control
use mpas_atmphys_driver
use mpas_atmphys_driver, only : physics_driver_before_dyn, physics_driver_after_dyn
use mpas_atmphys_manager
use mpas_atmphys_update
#endif
Expand All @@ -852,17 +852,39 @@ subroutine atm_do_timestep(domain, dt, itimestep)

call mpas_get_time(curr_time=currTime, dateTimeString=timeStamp, ierr=ierr)


!
! This is a call to all of the slow, column-oriented physics, except for resolved
! cloud microphysical processes. The microphysics is handled after the dynamics so
! that any supersaturation caused by various dynamical processes (advection, etc)
! may be removed before the end of each time step.
!
#ifdef DO_PHYSICS
!proceed with physics if moist_physics is set to true:
if(moist_physics) then
call physics_timetracker(domain,dt,clock,itimestep,xtime_s)
call physics_driver(domain,itimestep,xtime_s)
call physics_driver_before_dyn(domain,itimestep,xtime_s)
endif
#endif

call atm_timestep(domain, dt, timeStamp, itimestep)

!
! This part is the microphysics, which needs to happen after the dynamics
! time stepping.
!
#ifdef DO_PHYSICS
!proceed with physics if moist_physics is set to true:
if(moist_physics) then
call physics_driver_after_dyn(domain,itimestep,xtime_s)
endif
#endif

!
! This is the summary from the physics before the dynamics solver,
! the accumulated processing of the dynamics, and the final physics
! call after the dynamics completes.
call summarize_timestep(domain)

end subroutine atm_do_timestep


Expand Down
1 change: 1 addition & 0 deletions src/core_atmosphere/physics/Makefile
Original file line number Diff line number Diff line change
Expand Up @@ -79,6 +79,7 @@ mpas_atmphys_driver.o: \
mpas_atmphys_driver_radiation_sw.o \
mpas_atmphys_driver_sfclayer.o \
mpas_atmphys_driver_oml.o \
mpas_atmphys_driver_microphysics.o \
mpas_atmphys_constants.o \
mpas_atmphys_interface.o \
mpas_atmphys_update.o \
Expand Down
Loading