SCALE-RM
scale_urban_phy.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_stdio
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: urban_phy_setup
29 
30  abstract interface
31  subroutine urb( &
32  TR_URB_t, &
33  TB_URB_t, &
34  TG_URB_t, &
35  TC_URB_t, &
36  QC_URB_t, &
37  UC_URB_t, &
38  TRL_URB_t, &
39  TBL_URB_t, &
40  TGL_URB_t, &
41  RAINR_URB_t, &
42  RAINB_URB_t, &
43  RAING_URB_t, &
44  ROFF_URB_t, &
45  SFC_TEMP, &
46  ALBD_LW, &
47  ALBD_SW, &
48  MWFLX, &
49  MUFLX, &
50  MVFLX, &
51  SHFLX, &
52  LHFLX, &
53  GHFLX, &
54  Z0M, &
55  Z0H, &
56  Z0E, &
57  U10, &
58  V10, &
59  T2, &
60  Q2, &
61  TMPA, &
62  PRES, &
63  W1, &
64  U1, &
65  V1, &
66  DENS, &
67  QA, &
68  Z1, &
69  PBL, &
70  PRSS, &
71  LWD, &
72  SWD, &
73  PREC, &
74  TR_URB, &
75  TB_URB, &
76  TG_URB, &
77  TC_URB, &
78  QC_URB, &
79  UC_URB, &
80  TRL_URB, &
81  TBL_URB, &
82  TGL_URB, &
83  RAINR_URB, &
84  RAINB_URB, &
85  RAING_URB, &
86  ROFF_URB, &
87  LON, &
88  LAT, &
89  NOWDATE, &
90  dt )
91  use scale_precision
94  implicit none
95 
96  real(RP), intent(out) :: tr_urb_t (ia,ja)
97  real(RP), intent(out) :: tb_urb_t (ia,ja)
98  real(RP), intent(out) :: tg_urb_t (ia,ja)
99  real(RP), intent(out) :: tc_urb_t (ia,ja)
100  real(RP), intent(out) :: qc_urb_t (ia,ja)
101  real(RP), intent(out) :: uc_urb_t (ia,ja)
102  real(RP), intent(out) :: trl_urb_t (uks:uke,ia,ja)
103  real(RP), intent(out) :: tbl_urb_t (uks:uke,ia,ja)
104  real(RP), intent(out) :: tgl_urb_t (uks:uke,ia,ja)
105  real(RP), intent(out) :: rainr_urb_t(ia,ja)
106  real(RP), intent(out) :: rainb_urb_t(ia,ja)
107  real(RP), intent(out) :: raing_urb_t(ia,ja)
108  real(RP), intent(out) :: roff_urb_t (ia,ja)
109 
110  real(RP), intent(out) :: sfc_temp(ia,ja)
111  real(RP), intent(out) :: albd_lw (ia,ja)
112  real(RP), intent(out) :: albd_sw (ia,ja)
113  real(RP), intent(out) :: mwflx (ia,ja)
114  real(RP), intent(out) :: muflx (ia,ja)
115  real(RP), intent(out) :: mvflx (ia,ja)
116  real(RP), intent(out) :: shflx (ia,ja)
117  real(RP), intent(out) :: lhflx (ia,ja)
118  real(RP), intent(out) :: ghflx (ia,ja)
119  real(RP), intent(out) :: z0m (ia,ja)
120  real(RP), intent(out) :: z0h (ia,ja)
121  real(RP), intent(out) :: z0e (ia,ja)
122  real(RP), intent(out) :: u10 (ia,ja)
123  real(RP), intent(out) :: v10 (ia,ja)
124  real(RP), intent(out) :: t2 (ia,ja)
125  real(RP), intent(out) :: q2 (ia,ja)
126 
127  real(RP), intent(in) :: tmpa (ia,ja)
128  real(RP), intent(in) :: pres (ia,ja)
129  real(RP), intent(in) :: w1 (ia,ja)
130  real(RP), intent(in) :: u1 (ia,ja)
131  real(RP), intent(in) :: v1 (ia,ja)
132  real(RP), intent(in) :: dens (ia,ja)
133  real(RP), intent(in) :: qa (ia,ja)
134  real(RP), intent(in) :: z1 (ia,ja)
135  real(RP), intent(in) :: pbl (ia,ja)
136  real(RP), intent(in) :: prss (ia,ja)
137  real(RP), intent(in) :: lwd (ia,ja,2)
138  real(RP), intent(in) :: swd (ia,ja,2)
139  real(RP), intent(in) :: prec (ia,ja)
140 
141  real(RP), intent(in) :: tr_urb (ia,ja)
142  real(RP), intent(in) :: tb_urb (ia,ja)
143  real(RP), intent(in) :: tg_urb (ia,ja)
144  real(RP), intent(in) :: tc_urb (ia,ja)
145  real(RP), intent(in) :: qc_urb (ia,ja)
146  real(RP), intent(in) :: uc_urb (ia,ja)
147  real(RP), intent(in) :: trl_urb (uks:uke,ia,ja)
148  real(RP), intent(in) :: tbl_urb (uks:uke,ia,ja)
149  real(RP), intent(in) :: tgl_urb (uks:uke,ia,ja)
150  real(RP), intent(in) :: rainr_urb(ia,ja)
151  real(RP), intent(in) :: rainb_urb(ia,ja)
152  real(RP), intent(in) :: raing_urb(ia,ja)
153  real(RP), intent(in) :: roff_urb (ia,ja)
154 
155  real(RP), intent(in) :: lon
156  real(RP), intent(in) :: lat
157  integer, intent(in) :: nowdate(6)
158  real(DP), intent(in) :: dt
159  end subroutine urb
160  end interface
161  procedure(urb), pointer :: urban_phy => null()
162  public :: urban_phy
163 
164  !-----------------------------------------------------------------------------
165  !
166  !++ Public parameters & variables
167  !
168  !-----------------------------------------------------------------------------
169  !
170  !++ Private procedure
171  !
172  !-----------------------------------------------------------------------------
173  !
174  !++ Private parameters & variables
175  !
176  !-----------------------------------------------------------------------------
177 contains
178 
179  subroutine urban_phy_setup( URBAN_TYPE )
180  use scale_process, only: &
182  use scale_urban_phy_slc, only: &
185  implicit none
186 
187  character(len=*), intent(in) :: URBAN_TYPE
188  !---------------------------------------------------------------------------
189 
190  select case( urban_type )
191  case ( 'SLC' )
192  call urban_phy_slc_setup( urban_type )
194  case default
195  write(*,*) 'xxx invalid Urban type(', trim(urban_type), '). CHECK!'
196  call prc_mpistop
197  end select
198 
199  end subroutine urban_phy_setup
200 
201 end module scale_urban_phy
subroutine, public prc_mpistop
Abort MPI.
module URBAN / Physics
subroutine, public urban_phy_slc_setup(URBAN_TYPE)
Setup.
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
subroutine, public urban_phy_setup(URBAN_TYPE)
procedure(urb), pointer, public urban_phy
module grid index
integer, public ia
of x whole cells (local, with HALO)
module URBAN / Surface fluxes with Single-layer Canpoy Model
module PROCESS
subroutine, public urban_phy_slc(TR_URB_t, TB_URB_t, TG_URB_t, TC_URB_t, QC_URB_t, UC_URB_t, TRL_URB_t, TBL_URB_t, TGL_URB_t, RAINR_URB_t, RAINB_URB_t, RAING_URB_t, ROFF_URB_t, SFC_TEMP, ALBD_LW, ALBD_SW, MWFLX, MUFLX, MVFLX, SHFLX, LHFLX, GHFLX, Z0M, Z0H, Z0E, U10, V10, T2, Q2, TMPA, PRSA, W1, U1, V1, DENS, QA, Z1, PBL, PRSS, LWD, SWD, PREC, TR_URB, TB_URB, TG_URB, TC_URB, QC_URB, UC_URB, TRL_URB, TBL_URB, TGL_URB, RAINR_URB, RAINB_URB, RAING_URB, ROFF_URB, LON, LAT, NOWDATE, dt)
module PRECISION
module urban grid index
integer, public ja
of y whole cells (local, with HALO)