SCALE-RM
scale_atmos_phy_cp.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14  !------------------------------------------------------------------------------
15  !
16  !+++ used modules
17  !
18  use scale_precision
19  use scale_stdio
20  use scale_prof
22  use scale_tracer
23  private
24 
25  abstract interface
26  subroutine cp( &
27  DENS, &
28  MOMZ, &
29  MOMX, &
30  MOMY, &
31  RHOT, &
32  QTRC, &
33  DENS_t_CP, &
34  MOMZ_t_CP, &
35  MOMX_t_CP, &
36  MOMY_t_CP, &
37  RHOT_t_CP, &
38  RHOQ_t_CP, &
39  MFLX_cloudbase, &
40  SFLX_convrain, &
41  cloudtop, &
42  cloudbase, &
43  cldfrac_dp, &
44  cldfrac_sh, &
45  kf_nca, &
46  kf_w0avg )
47  use scale_precision
48  use scale_stdio
49  use scale_prof
51  use scale_tracer
52  real(RP), intent(in) :: dens(ka,ia,ja)
53  real(RP), intent(in) :: momx(ka,ia,ja)
54  real(RP), intent(in) :: momy(ka,ia,ja)
55  real(RP), intent(in) :: momz(ka,ia,ja)
56  real(RP), intent(in) :: rhot(ka,ia,ja)
57  real(RP), intent(in) :: qtrc(ka,ia,ja,qa)
58  real(RP), intent(inout) :: dens_t_cp(ka,ia,ja)
59  real(RP), intent(inout) :: momz_t_cp(ka,ia,ja)
60  real(RP), intent(inout) :: momx_t_cp(ka,ia,ja)
61  real(RP), intent(inout) :: momy_t_cp(ka,ia,ja)
62  real(RP), intent(inout) :: rhot_t_cp(ka,ia,ja)
63  real(RP), intent(inout) :: rhoq_t_cp(ka,ia,ja,qa)
64  real(RP), intent(inout) :: mflx_cloudbase(ia,ja)
65  real(RP), intent(inout) :: sflx_convrain(ia,ja)
66  real(RP), intent(inout) :: cloudtop(ia,ja)
67  real(RP), intent(inout) :: cloudbase(ia,ja)
68  real(RP), intent(inout) :: cldfrac_dp(ka,ia,ja)
69  real(RP), intent(inout) :: cldfrac_sh(ka,ia,ja)
70  real(RP), intent(inout) :: kf_nca(ia,ja)
71  real(RP), intent(inout) :: kf_w0avg(ka,ia,ja)
72  end subroutine cp
73  end interface
74  procedure(cp), pointer :: atmos_phy_cp => null()
75  public :: atmos_phy_cp
76  public :: atmos_phy_cp_setup
77 
78  !-----------------------------------------------------------------------------
79 contains
80  !-----------------------------------------------------------------------------
82  !------------------------------------------------------------------------------
83  subroutine atmos_phy_cp_setup( CP_TYPE )
84  use scale_process, only: &
86 #define EXTM(pre, name, post) pre ## name ## post
87 #define NAME(pre, name, post) EXTM(pre, name, post)
88 #ifdef CP
89  use name(scale_atmos_phy_mp_, cp,), only: &
90  name(atmos_phy_cp_, cp, _setup), &
91  name(atmos_phy_cp_, cp,), &
92 #else
93  use scale_atmos_phy_cp_kf, only: &
96 #endif
97  implicit none
98 
99  character(len=*), intent(in) :: CP_TYPE
100  !------------------------------------------------------------------------------
101 
102 #ifdef CP
103  call name(atmos_phy_cp_, cp, _setup)( cp_type )
104  atmos_phy_cp => name(atmos_phy_cp_, cp,)
105 #else
106  select case ( cp_type )
107  case( 'KF' )
108  call atmos_phy_cp_kf_setup( cp_type )
110  end select
111 #endif
112 
113  return
114  end subroutine atmos_phy_cp_setup
115 
116 end module scale_atmos_phy_cp
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_phy_cp_setup(CP_TYPE)
Setup Cumulus parameterization.
module ATMOSPHERE / Physics Cumulus Parameterization
subroutine, public atmos_phy_cp_kf_setup(CP_TYPE)
Setup.
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
module PROCESS
module ATMOSPHERE / Physics Cumulus Parameterization
procedure(cp), pointer, public atmos_phy_cp
module profiler
Definition: scale_prof.F90:10
module PRECISION
subroutine, public atmos_phy_cp_kf(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, DENS_t_CP, MOMZ_t_CP, MOMX_t_CP, MOMY_t_CP, RHOT_t_CP, RHOQ_t_CP, MFLX_cloudbase, SFLX_convrain, cloudtop, cloudbase, cldfrac_dp, cldfrac_sh, nca, w0avg)
integer, public ja
of y whole cells (local, with HALO)