SCALE-RM
scale_atmos_phy_tb.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
16  !-----------------------------------------------------------------------------
17  !
18  !++ used modules
19  !
20  use scale_precision
21  use scale_stdio
22  use scale_prof
24  use scale_tracer
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public procedure
31  !
32  public :: atmos_phy_tb_setup
33 
34  !-----------------------------------------------------------------------------
35  !
36  !++ Public parameters & variables
37  !
38  !-----------------------------------------------------------------------------
39  !
40  !++ Private procedure
41  !
42  !-----------------------------------------------------------------------------
43  !
44  !++ Private parameters & variables
45  !
46  abstract interface
47  subroutine tb( &
48  qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, & ! (out)
49  qflx_sgs_rhot, qflx_sgs_rhoq, & ! (out)
50  tke, & ! (inout)
51  tke_t, nu_c, ri, pr, n2, & ! (out) diagnostic variables
52  momz, momx, momy, rhot, dens, qtrc, & ! (in)
53  sflx_mw, sflx_mu, sflx_mv, sflx_sh, sflx_qv, & ! (in)
54  gsqrt, j13g, j23g, j33g, mapf, dt ) ! (in)
55  use scale_precision
57  use scale_tracer
58  implicit none
59  real(RP), intent(out) :: qflx_sgs_momz(ka,ia,ja,3)
60  real(RP), intent(out) :: qflx_sgs_momx(ka,ia,ja,3)
61  real(RP), intent(out) :: qflx_sgs_momy(ka,ia,ja,3)
62  real(RP), intent(out) :: qflx_sgs_rhot(ka,ia,ja,3)
63  real(RP), intent(out) :: qflx_sgs_rhoq(ka,ia,ja,3,qa)
64 
65  real(RP), intent(inout) :: tke(ka,ia,ja)
66  real(RP), intent(out) :: tke_t(ka,ia,ja) ! tendency TKE
67  real(RP), intent(out) :: nu_c(ka,ia,ja) ! eddy viscosity (center)
68  real(RP), intent(out) :: ri (ka,ia,ja) ! Richardson number
69  real(RP), intent(out) :: pr (ka,ia,ja) ! Prantle number
70  real(RP), intent(out) :: n2 (ka,ia,ja) ! squared Brunt-Vaisala frequency
71 
72  real(RP), intent(in) :: momz(ka,ia,ja)
73  real(RP), intent(in) :: momx(ka,ia,ja)
74  real(RP), intent(in) :: momy(ka,ia,ja)
75  real(RP), intent(in) :: rhot(ka,ia,ja)
76  real(RP), intent(in) :: dens(ka,ia,ja)
77  real(RP), intent(in) :: qtrc(ka,ia,ja,qa)
78 
79  real(RP), intent(in) :: sflx_mw(ia,ja)
80  real(RP), intent(in) :: sflx_mu(ia,ja)
81  real(RP), intent(in) :: sflx_mv(ia,ja)
82  real(RP), intent(in) :: sflx_sh(ia,ja)
83  real(RP), intent(in) :: sflx_qv(ia,ja)
84 
85  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
86  real(RP), intent(in) :: j13g (ka,ia,ja,7)
87  real(RP), intent(in) :: j23g (ka,ia,ja,7)
88  real(RP), intent(in) :: j33g
89  real(RP), intent(in) :: mapf (ia,ja,2,4)
90  real(DP), intent(in) :: dt
91  end subroutine tb
92  end interface
93  procedure(tb), pointer :: atmos_phy_tb => null()
94  public :: atmos_phy_tb
95 
96  !-----------------------------------------------------------------------------
97 contains
98 
99  subroutine atmos_phy_tb_setup( &
100  TB_TYPE, &
101  CDZ, CDX, CDY, &
102  CZ )
103  use scale_process, only: &
105 #define EXTM(pre, name, post) pre ## name ## post
106 #define NAME(pre, name, post) EXTM(pre, name, post)
107 #ifdef TB
108  use name(scale_atmos_phy_tb_, tb,), only: &
109  name(atmos_phy_tb_, tb, _setup), &
110  name(atmos_phy_tb_, tb,)
111 #else
112  use scale_atmos_phy_tb_smg, only: &
115  use scale_atmos_phy_tb_d1980, only: &
118  use scale_atmos_phy_tb_dns, only: &
121  use scale_atmos_phy_tb_mynn, only: &
124  use scale_atmos_phy_tb_hybrid, only: &
127  use scale_atmos_phy_tb_dummy, only: &
130 #endif
131  implicit none
132 
133  character(len=*), intent(in) :: TB_TYPE
134 
135  real(RP), intent(in) :: CDZ(ka)
136  real(RP), intent(in) :: CDX(ia)
137  real(RP), intent(in) :: CDY(ja)
138  real(RP), intent(in) :: CZ (ka,ia,ja)
139  !---------------------------------------------------------------------------
140 
141 #ifdef TB
142  call name(atmos_phy_tb_, tb, _setup)( &
143  tb_type, &
144  cdz, cdx, cdy, &
145  cz )
146  atmos_phy_tb => name(atmos_phy_tb_, tb, )
147 #else
148  select case( tb_type )
149  case ( 'SMAGORINSKY' )
150  call atmos_phy_tb_smg_setup( &
151  tb_type, &
152  cdz, cdx, cdy, &
153  cz )
155 
156  case ( 'D1980' )
158  tb_type, &
159  cdz, cdx, cdy, &
160  cz )
162 
163  case ( 'DNS' )
164  call atmos_phy_tb_dns_setup( &
165  tb_type, &
166  cdz, cdx, cdy, &
167  cz )
169 
170  case ( 'MYNN' )
172  tb_type, &
173  cdz, cdx, cdy, &
174  cz )
176 
177  case ('HYBRID')
179  tb_type, &
180  cdz, cdx, cdy, &
181  cz )
183 
184  case ('OFF')
185 
186  ! do nothing
187 
188  case default
189 
190  write(*,*) 'xxx ATMOS_PHY_TB_TYPE is invalid'
191  call prc_mpistop
192 
193  end select
194 #endif
195 
196  return
197  end subroutine atmos_phy_tb_setup
198 
199 end module scale_atmos_phy_tb
module ATMOSPHERE / Physics Turbulence
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_phy_tb_d1980_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
procedure(tb), pointer, public atmos_phy_tb
subroutine, public atmos_phy_tb_mynn(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, Nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
subroutine, public atmos_phy_tb_d1980(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, Km, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
subroutine, public atmos_phy_tb_dns(qflx_sgs_MOMZ, qflx_sgs_MOMX, qflx_sgs_MOMY, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
subroutine, public atmos_phy_tb_mynn_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
module grid index
module TRACER
module ATMOSPHERE / Physics Turbulence
integer, public ia
of x whole cells (local, with HALO)
module ATMOSPHERE / Physics Turbulence
integer, public ka
of z whole cells (local, with HALO)
subroutine, public atmos_phy_tb_dns_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
subroutine, public atmos_phy_tb_smg(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
subroutine, public atmos_phy_tb_smg_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
subroutine, public atmos_phy_tb_dummy_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
module ATMOSPHERE / Physics Turbulence
module PROCESS
subroutine, public atmos_phy_tb_setup(TB_TYPE, CDZ, CDX, CDY, CZ)
module ATMOSPHERE / Physics Turbulence
subroutine, public atmos_phy_tb_hybrid_setup(TB_TYPE, CDZ, CDX, CDY, CZ)
module profiler
Definition: scale_prof.F90:10
module ATMOSPHERE / Physics Turbulence
module PRECISION
subroutine, public atmos_phy_tb_hybrid(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, Nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
module ATMOSPHERE / Physics Turbulence
subroutine, public atmos_phy_tb_dummy(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, nu_C, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, sflx_mw, sflx_mu, sflx_mv, sflx_sh, sflx_qv, GSQRT, J13G, J23G, J33G, MAPF, dt)
integer, public ja
of y whole cells (local, with HALO)