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  abstract interface
33  subroutine tb( &
34  qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, &
35  qflx_sgs_rhot, qflx_sgs_rhoq, &
36  RHOQ_t, &
37  nu_C, Ri, Pr, &
38  MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, &
39  SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, &
40  GSQRT, J13G, J23G, J33G, MAPF, dt )
41  use scale_precision
43  use scale_tracer
44  implicit none
45 
46  real(RP), intent(out) :: qflx_sgs_momz(KA,IA,JA,3)
47  real(RP), intent(out) :: qflx_sgs_momx(KA,IA,JA,3)
48  real(RP), intent(out) :: qflx_sgs_momy(KA,IA,JA,3)
49  real(RP), intent(out) :: qflx_sgs_rhot(KA,IA,JA,3)
50  real(RP), intent(out) :: qflx_sgs_rhoq(KA,IA,JA,3,QA)
51 
52  real(RP), intent(inout) :: RHOQ_t (KA,IA,JA,QA) ! tendency of rho * QTRC
53 
54  real(RP), intent(out) :: nu_C (KA,IA,JA) ! eddy viscosity (center)
55  real(RP), intent(out) :: Ri (KA,IA,JA) ! Richardson number
56  real(RP), intent(out) :: Pr (KA,IA,JA) ! Prantle number
57 
58  real(RP), intent(in) :: MOMZ (KA,IA,JA)
59  real(RP), intent(in) :: MOMX (KA,IA,JA)
60  real(RP), intent(in) :: MOMY (KA,IA,JA)
61  real(RP), intent(in) :: RHOT (KA,IA,JA)
62  real(RP), intent(in) :: DENS (KA,IA,JA)
63  real(RP), intent(in) :: QTRC (KA,IA,JA,QA)
64  real(RP), intent(in) :: N2 (KA,IA,JA)
65 
66  real(RP), intent(in) :: SFLX_MW (IA,JA)
67  real(RP), intent(in) :: SFLX_MU (IA,JA)
68  real(RP), intent(in) :: SFLX_MV (IA,JA)
69  real(RP), intent(in) :: SFLX_SH (IA,JA)
70  real(RP), intent(in) :: SFLX_Q (IA,JA,QA)
71 
72  real(RP), intent(in) :: GSQRT (KA,IA,JA,7)
73  real(RP), intent(in) :: J13G (KA,IA,JA,7)
74  real(RP), intent(in) :: J23G (KA,IA,JA,7)
75  real(RP), intent(in) :: J33G
76  real(RP), intent(in) :: MAPF (IA,JA,2,4)
77  real(DP), intent(in) :: dt
78  end subroutine tb
79 
80  subroutine su( &
81  CDZ, CDX, CDY, CZ )
82  use scale_precision
84  use scale_tracer
85  implicit none
86 
87  real(RP), intent(in) :: CDZ(KA)
88  real(RP), intent(in) :: CDX(IA)
89  real(RP), intent(in) :: CDY(JA)
90  real(RP), intent(in) :: CZ (KA,IA,JA)
91  end subroutine su
92  end interface
93 
94  procedure(tb), pointer :: atmos_phy_tb => null()
95 
96  procedure(su), pointer :: atmos_phy_tb_setup => null()
97 
98  public :: atmos_phy_tb_config
99  public :: atmos_phy_tb
100  public :: atmos_phy_tb_setup
101 
102  !-----------------------------------------------------------------------------
103  !
104  !++ Public parameters & variables
105  !
106  integer, public :: i_tke = -1
107 
108  !-----------------------------------------------------------------------------
109  !
110  !++ Private procedure
111  !
112  !-----------------------------------------------------------------------------
113  !
114  !++ Private parameters & variables
115  !
116  !-----------------------------------------------------------------------------
117 contains
118  !-----------------------------------------------------------------------------
120  subroutine atmos_phy_tb_config( TB_TYPE )
121  use scale_process, only: &
123  use scale_atmos_phy_tb_smg, only: &
127  use scale_atmos_phy_tb_d1980, only: &
131  use scale_atmos_phy_tb_dns, only: &
135  use scale_atmos_phy_tb_mynn, only: &
139  use scale_atmos_phy_tb_hybrid, only: &
143  implicit none
144 
145  character(len=*), intent(in) :: tb_type
146  !---------------------------------------------------------------------------
147 
148  if( io_l ) write(io_fid_log,*) '*** => ', trim(tb_type), ' is selected.'
149 
150  select case( tb_type )
151  case( 'SMAGORINSKY' )
152 
153  call atmos_phy_tb_smg_config( tb_type, & ! [IN]
154  i_tke ) ! [OUT]
155 
158 
159  case( 'D1980' )
160 
161  call atmos_phy_tb_d1980_config( tb_type, & ! [IN]
162  i_tke ) ! [OUT]
163 
166 
167  case( 'DNS' )
168 
169  call atmos_phy_tb_dns_config( tb_type, & ! [IN]
170  i_tke ) ! [OUT]
171 
174 
175  case( 'MYNN' )
176 
177  call atmos_phy_tb_mynn_config( tb_type, & ! [IN]
178  i_tke ) ! [OUT]
179 
182 
183  case('HYBRID')
184 
185  call atmos_phy_tb_hybrid_config( tb_type, & ! [IN]
186  i_tke ) ! [OUT]
187 
190 
191  case('OFF')
192 
193  ! do nothing
194 
195  case default
196 
197  write(*,*) 'xxx ATMOS_PHY_TB_TYPE is invalid: ', tb_type
198  call prc_mpistop
199 
200  end select
201 
202  return
203  end subroutine atmos_phy_tb_config
204 
205 end module scale_atmos_phy_tb
module ATMOSPHERE / Physics Turbulence
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_phy_tb_d1980(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, RHOQ_t, Km, Ri, Pr, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, 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, RHOQ_t, nu, Ri, Pr, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, GSQRT, J13G, J23G, J33G, MAPF, dt)
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
subroutine, public atmos_phy_tb_d1980_config(TYPE_TB, I_TKE_out)
Config.
procedure(tb), pointer, public atmos_phy_tb
module STDIO
Definition: scale_stdio.F90:12
subroutine, public atmos_phy_tb_smg_setup(CDZ, CDX, CDY, CZ)
Setup.
subroutine, public atmos_phy_tb_mynn_setup(CDZ, CDX, CDY, CZ)
Setup.
module grid index
module TRACER
module ATMOSPHERE / Physics Turbulence
subroutine, public atmos_phy_tb_dns_setup(CDZ, CDX, CDY, CZ)
module ATMOSPHERE / Physics Turbulence
subroutine, public atmos_phy_tb_hybrid_config(TB_TYPE, I_TKE_out)
Config.
subroutine, public atmos_phy_tb_mynn(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, RHOQ_t, Nu, Ri, Pr, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2_in, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, GSQRT, J13G, J23G, J33G, MAPF, dt)
subroutine, public atmos_phy_tb_dns_config(TYPE_TB, I_TKE_out)
Config.
subroutine, public atmos_phy_tb_mynn_config(TYPE_TB, I_TKE_out)
Config.
module ATMOSPHERE / Physics Turbulence
module PROCESS
subroutine, public atmos_phy_tb_hybrid_setup(CDZ, CDX, CDY, CZ)
Setup.
module ATMOSPHERE / Physics Turbulence
subroutine, public atmos_phy_tb_smg(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, RHOQ_t, Nu, Ri, Pr, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, GSQRT, J13G, J23G, J33G, MAPF, dt)
module profiler
Definition: scale_prof.F90:10
subroutine, public atmos_phy_tb_d1980_setup(CDZ, CDX, CDY, CZ)
Setup.
module PRECISION
subroutine, public atmos_phy_tb_hybrid(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, RHOQ_t, Nu, Ri, Pr, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, GSQRT, J13G, J23G, J33G, MAPF, dt)
procedure(su), pointer, public atmos_phy_tb_setup
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public atmos_phy_tb_config(TB_TYPE)
register
module ATMOSPHERE / Physics Turbulence
subroutine, public atmos_phy_tb_smg_config(TYPE_TB, I_TKE_out)
Config.