SCALE-RM
Functions/Subroutines
mod_atmos_phy_sf_driver Module Reference

module ATMOSPHERE / Physics Surface fluxes More...

Functions/Subroutines

subroutine, public atmos_phy_sf_driver_setup
 Setup. More...
 
subroutine, public atmos_phy_sf_driver_resume
 Resume. More...
 
subroutine, public atmos_phy_sf_driver (update_flag)
 Driver. More...
 

Detailed Description

module ATMOSPHERE / Physics Surface fluxes

Description
Flux from/to bottom boundary of atmosphere (surface)
Author
Team SCALE
History
  • 2013-12-05 (S.Nishizawa) [new]
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
GHFLX ground heat flux (merged) W/m2 SFLX_GH
LHFLX latent heat flux (merged) W/m2 SFLX_LH
MSLP mean sea-level pressure Pa MSLP
MUFLX u-momentum flux (merged) kg/m/s2 SFLX_MU
MVFLX v-momentum flux (merged) kg/m/s2 SFLX_MV
MWFLX w-momentum flux (merged) kg/m/s2 SFLX_MW
Q2 2m specific humidity kg/kg Q2
SFC_ALB_LW 'surface albedo (longwave merged)' SFC_albedo
SFC_ALB_SW 'surface albedo (shortwave merged)' SFC_albedo
SFC_DENS surface atmospheric density kg/m3 SFC_DENS
SFC_PRES surface atmospheric pressure Pa SFC_PRES
SFC_TEMP surface skin temperature (merged) K SFC_TEMP
SFC_Z0E roughness length (vapor) m SFC_Z0E
SFC_Z0H roughness length (heat) m SFC_Z0H
SFC_Z0M roughness length (momentum) m SFC_Z0M
SHFLX sensible heat flux (merged) W/m2 SFLX_SH
T2 2m air temperature K T2
U10 10m x-wind m/s U10
Uabs10 10m absolute wind m/s Uabs10
V10 10m y-wind m/s V10

Function/Subroutine Documentation

◆ atmos_phy_sf_driver_setup()

subroutine, public mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup ( )

Setup.

Definition at line 56 of file mod_atmos_phy_sf_driver.f90.

References scale_atmos_phy_sf::atmos_phy_sf_setup(), mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0e, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0h, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0m, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_lh, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mw, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, mod_atmos_admin::atmos_phy_sf_type, mod_atmos_admin::atmos_sw_phy_sf, mod_cpl_admin::cpl_sw, scale_stdio::io_fid_log, scale_stdio::io_l, and scale_process::prc_mpistop().

Referenced by mod_atmos_driver::atmos_driver_setup().

56  use scale_process, only: &
58  use scale_atmos_phy_sf, only: &
60  use mod_atmos_admin, only: &
63  use mod_atmos_phy_sf_vars, only: &
64  sfc_z0m => atmos_phy_sf_sfc_z0m, &
65  sfc_z0h => atmos_phy_sf_sfc_z0h, &
66  sfc_z0e => atmos_phy_sf_sfc_z0e, &
67  sflx_mw => atmos_phy_sf_sflx_mw, &
68  sflx_mu => atmos_phy_sf_sflx_mu, &
69  sflx_mv => atmos_phy_sf_sflx_mv, &
70  sflx_sh => atmos_phy_sf_sflx_sh, &
71  sflx_lh => atmos_phy_sf_sflx_lh, &
72  sflx_qtrc => atmos_phy_sf_sflx_qtrc
73  use mod_cpl_admin, only: &
74  cpl_sw
75  implicit none
76  !---------------------------------------------------------------------------
77 
78  if( io_l ) write(io_fid_log,*)
79  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[ATMOS PHY_SF] / Origin[SCALE-RM]'
80 
81  if ( atmos_sw_phy_sf ) then
82 
83  ! setup library component
85 
86  if ( .NOT. cpl_sw ) then
87  if( io_l ) write(io_fid_log,*) '*** Coupler is disabled.'
88  endif
89 
90  else
91 
92  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
93  if( io_l ) write(io_fid_log,*) '*** surface fluxes are set to zero.'
94  sflx_mw(:,:) = 0.0_rp
95  sflx_mu(:,:) = 0.0_rp
96  sflx_mv(:,:) = 0.0_rp
97  sflx_sh(:,:) = 0.0_rp
98  sflx_lh(:,:) = 0.0_rp
99  sflx_qtrc(:,:,:) = 0.0_rp
100 
101  if( io_l ) write(io_fid_log,*) '*** SFC_TEMP, SFC_albedo is set in ATMOS_PHY_SF_vars.'
102 
103  endif
104 
105  return
module ATMOS admin
subroutine, public atmos_phy_sf_setup(SF_TYPE)
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sflx_qtrc
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
character(len=h_short), public atmos_phy_sf_type
module ATMOSPHERIC Surface Variables
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
module PROCESS
logical, public atmos_sw_phy_sf
module ATMOSPHERE / Physics Surface fluxes
logical, public cpl_sw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
module Coupler admin
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_driver_resume()

subroutine, public mod_atmos_phy_sf_driver::atmos_phy_sf_driver_resume ( )

Resume.

Definition at line 111 of file mod_atmos_phy_sf_driver.f90.

References atmos_phy_sf_driver(), mod_atmos_admin::atmos_sw_phy_sf, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by mod_atmos_driver::atmos_driver_resume2().

111  use mod_atmos_admin, only: &
113  implicit none
114 
115  if ( atmos_sw_phy_sf ) then
116 
117  ! run once (only for the diagnostic value)
118  call prof_rapstart('ATM_SurfaceFlux', 1)
119  call atmos_phy_sf_driver( update_flag = .true. )
120  call prof_rapend ('ATM_SurfaceFlux', 1)
121 
122  end if
123 
124  return
module ATMOS admin
logical, public atmos_sw_phy_sf
subroutine, public atmos_phy_sf_driver(update_flag)
Driver.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_driver()

subroutine, public mod_atmos_phy_sf_driver::atmos_phy_sf_driver ( logical, intent(in)  update_flag)

Driver.

Definition at line 130 of file mod_atmos_phy_sf_driver.f90.

References scale_atmos_bottom::atmos_bottom_estimate(), mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_lw_dn, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_sw_dn, scale_atmos_phy_sf::atmos_phy_sf, mod_atmos_phy_sf_vars::atmos_phy_sf_dens_t, mod_atmos_phy_sf_vars::atmos_phy_sf_momx_t, mod_atmos_phy_sf_vars::atmos_phy_sf_momy_t, mod_atmos_phy_sf_vars::atmos_phy_sf_momz_t, mod_atmos_phy_sf_vars::atmos_phy_sf_q2, mod_atmos_phy_sf_vars::atmos_phy_sf_rhoq_t, mod_atmos_phy_sf_vars::atmos_phy_sf_rhot_t, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_albedo, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_dens, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_pres, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_temp, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0e, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0h, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0m, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_gh, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_lh, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mw, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, mod_atmos_phy_sf_vars::atmos_phy_sf_t2, mod_atmos_phy_sf_vars::atmos_phy_sf_u10, mod_atmos_phy_sf_vars::atmos_phy_sf_v10, mod_cpl_admin::cpl_sw, mod_atmos_vars::dens, mod_atmos_vars::dens_av, mod_atmos_vars::dens_tp, scale_grid::grid_rcdz, scale_grid::grid_rfdz, scale_gridtrans::gtrans_gsqrt, scale_atmos_hydrometeor::i_qv, scale_gridtrans::i_uyz, scale_gridtrans::i_xvz, scale_gridtrans::i_xyw, scale_gridtrans::i_xyz, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ks, mod_atmos_vars::momx_tp, mod_atmos_vars::momy_tp, mod_atmos_vars::momz_tp, mod_atmos_vars::pres, scale_tracer::qa, mod_atmos_vars::qtrc, mod_atmos_vars::qtrc_av, scale_grid_real::real_cz, scale_grid_real::real_z1, mod_atmos_vars::rhoq_tp, mod_atmos_vars::rhot, mod_atmos_vars::rhot_av, mod_atmos_vars::rhot_tp, scale_rm_statistics::statistics_checktotal, mod_atmos_vars::temp, scale_time::time_dtsec_atmos_phy_sf, scale_topography::topo_zsfc, scale_tracer::tracer_cp, scale_tracer::tracer_mass, scale_tracer::tracer_name, scale_tracer::tracer_r, mod_atmos_vars::u, mod_atmos_vars::v, and mod_atmos_vars::w.

Referenced by mod_atmos_driver::atmos_driver(), and atmos_phy_sf_driver_resume().

130  use scale_comm, only: &
131  comm_vars8, &
132  comm_wait
133  use scale_grid, only: &
134  rcdz => grid_rcdz, &
135  rfdz => grid_rfdz
136  use scale_grid_real, only: &
137  real_cz, &
138  real_z1
139  use scale_gridtrans, only: &
140  gsqrt => gtrans_gsqrt, &
141  i_xyz, &
142  i_uyz, &
143  i_xvz, &
144  i_xyw
145  use scale_topography, only: &
146  topo_zsfc
147  use scale_time, only: &
148  dt_sf => time_dtsec_atmos_phy_sf
149  use scale_rm_statistics, only: &
151  stat_total
152  use scale_history, only: &
153  hist_in
154  use scale_atmos_bottom, only: &
155  bottom_estimate => atmos_bottom_estimate
156  use scale_atmos_hydrostatic, only: &
157  barometric_law_mslp => atmos_hydrostatic_barometric_law_mslp
158  use scale_atmos_thermodyn, only: &
159  thermodyn_qd => atmos_thermodyn_qd, &
160  thermodyn_cp => atmos_thermodyn_cp, &
161  thermodyn_r => atmos_thermodyn_r
162  use scale_atmos_hydrometeor, only: &
163  i_qv
164  use scale_atmos_phy_sf, only: &
166  use mod_atmos_vars, only: &
167  dens => dens_av, &
168  rhot => rhot_av, &
169  qtrc => qtrc_av, &
170  temp, &
171  pres, &
172  w, &
173  u, &
174  v, &
175  dens_t => dens_tp, &
176  momz_t => momz_tp, &
177  momx_t => momx_tp, &
178  momy_t => momy_tp, &
179  rhot_t => rhot_tp, &
180  rhoq_t => rhoq_tp
181  use mod_atmos_phy_rd_vars, only: &
182  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
183  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn
184  use mod_atmos_phy_sf_vars, only: &
185  dens_t_sf => atmos_phy_sf_dens_t, &
186  momz_t_sf => atmos_phy_sf_momz_t, &
187  momx_t_sf => atmos_phy_sf_momx_t, &
188  momy_t_sf => atmos_phy_sf_momy_t, &
189  rhot_t_sf => atmos_phy_sf_rhot_t, &
190  rhoq_t_sf => atmos_phy_sf_rhoq_t, &
191  sfc_dens => atmos_phy_sf_sfc_dens, &
192  sfc_pres => atmos_phy_sf_sfc_pres, &
193  sfc_temp => atmos_phy_sf_sfc_temp, &
194  sfc_albedo => atmos_phy_sf_sfc_albedo, &
195  sfc_z0m => atmos_phy_sf_sfc_z0m, &
196  sfc_z0h => atmos_phy_sf_sfc_z0h, &
197  sfc_z0e => atmos_phy_sf_sfc_z0e, &
198  sflx_mw => atmos_phy_sf_sflx_mw, &
199  sflx_mu => atmos_phy_sf_sflx_mu, &
200  sflx_mv => atmos_phy_sf_sflx_mv, &
201  sflx_sh => atmos_phy_sf_sflx_sh, &
202  sflx_lh => atmos_phy_sf_sflx_lh, &
203  sflx_gh => atmos_phy_sf_sflx_gh, &
204  sflx_qtrc => atmos_phy_sf_sflx_qtrc, &
205  u10 => atmos_phy_sf_u10, &
206  v10 => atmos_phy_sf_v10, &
207  t2 => atmos_phy_sf_t2, &
208  q2 => atmos_phy_sf_q2
209  use mod_cpl_admin, only: &
210  cpl_sw
211  implicit none
212 
213  logical, intent(in) :: update_flag
214 
215  real(RP) :: Uabs10(IA,JA) ! 10m absolute wind [m/s]
216  real(RP) :: MSLP (IA,JA) ! mean sea-level pressure [Pa]
217  real(RP) :: total ! dummy
218 
219  real(RP) :: q(QA)
220  real(RP) :: qdry
221  real(RP) :: Rtot
222  real(RP) :: CPtot
223 
224  real(RP) :: work
225 
226  integer :: i, j, iq
227  !---------------------------------------------------------------------------
228 
229  if ( update_flag ) then
230 
231  ! update surface density, surface pressure
232  call bottom_estimate( dens(:,:,:), & ! [IN]
233  pres(:,:,:), & ! [IN]
234  real_cz(:,:,:), & ! [IN]
235  topo_zsfc(:,:), & ! [IN]
236  real_z1(:,:), & ! [IN]
237  sfc_dens(:,:), & ! [OUT]
238  sfc_pres(:,:) ) ! [OUT]
239 
240  if ( .NOT. cpl_sw ) then
241  call atmos_phy_sf( temp(ks,:,:), & ! [IN]
242  pres(ks,:,:), & ! [IN]
243  w(ks,:,:), & ! [IN]
244  u(ks,:,:), & ! [IN]
245  v(ks,:,:), & ! [IN]
246  dens(ks,:,:), & ! [IN]
247  qtrc(ks,:,:,:), & ! [IN]
248  real_z1(:,:), & ! [IN]
249  dt_sf, & ! [IN]
250  sfc_dens(:,:), & ! [IN]
251  sfc_pres(:,:), & ! [IN]
252  sflx_lw_dn(:,:), & ! [IN]
253  sflx_sw_dn(:,:), & ! [IN]
254  sfc_temp(:,:), & ! [IN]
255  sfc_albedo(:,:,:), & ! [IN]
256  sfc_z0m(:,:), & ! [INOUT]
257  sfc_z0h(:,:), & ! [INOUT]
258  sfc_z0e(:,:), & ! [INOUT]
259  sflx_mw(:,:), & ! [OUT]
260  sflx_mu(:,:), & ! [OUT]
261  sflx_mv(:,:), & ! [OUT]
262  sflx_sh(:,:), & ! [OUT]
263  sflx_lh(:,:), & ! [OUT]
264  sflx_qtrc(:,:,:), & ! [OUT]
265  u10(:,:), & ! [OUT]
266  v10(:,:), & ! [OUT]
267  t2(:,:), & ! [OUT]
268  q2(:,:) ) ! [OUT]
269  endif
270 
271 !OCL XFILL
272  do j = js, je
273  do i = is, ie
274  uabs10(i,j) = sqrt( u10(i,j)**2 + v10(i,j)**2 )
275  end do
276  end do
277 
278  call barometric_law_mslp( mslp(:,:), sfc_pres(:,:), t2(:,:), topo_zsfc(:,:) )
279 
280  call hist_in( sfc_dens(:,:), 'SFC_DENS', 'surface atmospheric density', 'kg/m3' )
281  call hist_in( sfc_pres(:,:), 'SFC_PRES', 'surface atmospheric pressure', 'Pa' )
282  call hist_in( sfc_temp(:,:), 'SFC_TEMP', 'surface skin temperature (merged)', 'K' )
283  call hist_in( sfc_albedo(:,:,i_lw), 'SFC_ALB_LW', 'surface albedo (longwave,merged)', '1' , nohalo=.true. )
284  call hist_in( sfc_albedo(:,:,i_sw), 'SFC_ALB_SW', 'surface albedo (shortwave,merged)', '1' , nohalo=.true. )
285  call hist_in( sfc_z0m(:,:), 'SFC_Z0M', 'roughness length (momentum)', 'm' , nohalo=.true. )
286  call hist_in( sfc_z0h(:,:), 'SFC_Z0H', 'roughness length (heat)', 'm' , nohalo=.true. )
287  call hist_in( sfc_z0e(:,:), 'SFC_Z0E', 'roughness length (vapor)', 'm' , nohalo=.true. )
288  call hist_in( sflx_mw(:,:), 'MWFLX', 'w-momentum flux (merged)', 'kg/m/s2' )
289  call hist_in( sflx_mu(:,:), 'MUFLX', 'u-momentum flux (merged)', 'kg/m/s2' )
290  call hist_in( sflx_mv(:,:), 'MVFLX', 'v-momentum flux (merged)', 'kg/m/s2' )
291  call hist_in( sflx_sh(:,:), 'SHFLX', 'sensible heat flux (merged)', 'W/m2' , nohalo=.true. )
292  call hist_in( sflx_lh(:,:), 'LHFLX', 'latent heat flux (merged)', 'W/m2' , nohalo=.true. )
293  call hist_in( sflx_gh(:,:), 'GHFLX', 'ground heat flux (merged)', 'W/m2' , nohalo=.true. )
294  call hist_in( uabs10(:,:), 'Uabs10', '10m absolute wind', 'm/s' , nohalo=.true. )
295  call hist_in( u10(:,:), 'U10', '10m x-wind', 'm/s' , nohalo=.true. )
296  call hist_in( v10(:,:), 'V10', '10m y-wind', 'm/s' , nohalo=.true. )
297  call hist_in( t2(:,:), 'T2 ', '2m air temperature', 'K' , nohalo=.true. )
298  call hist_in( q2(:,:), 'Q2 ', '2m specific humidity', 'kg/kg' , nohalo=.true. )
299  call hist_in( mslp(:,:), 'MSLP', 'mean sea-level pressure', 'Pa' )
300 
301  call comm_vars8( sflx_mu(:,:), 1 )
302  call comm_vars8( sflx_mv(:,:), 2 )
303  call comm_wait ( sflx_mu(:,:), 1 )
304  call comm_wait ( sflx_mv(:,:), 2 )
305 
306 !OCL XFILL
307  do j = js, je
308  do i = is, ie
309  momz_t_sf(i,j) = sflx_mw(i,j) * rfdz(ks) / gsqrt(ks,i,j,i_xyw)
310  enddo
311  enddo
312 
313 !OCL XFILL
314  do j = js, je
315  do i = is, ie
316  momx_t_sf(i,j) = 0.5_rp * ( sflx_mu(i,j) + sflx_mu(i+1,j) ) * rcdz(ks) / gsqrt(ks,i,j,i_uyz)
317  enddo
318  enddo
319 
320 !OCL XFILL
321  do j = js, je
322  do i = is, ie
323  momy_t_sf(i,j) = 0.5_rp * ( sflx_mv(i,j) + sflx_mv(i,j+1) ) * rcdz(ks) / gsqrt(ks,i,j,i_xvz)
324  enddo
325  enddo
326 
327  do j = js, je
328  do i = is, ie
329  do iq = 1, qa
330  q(iq) = qtrc(ks,i,j,iq)
331  enddo
332  call thermodyn_qd( qdry, q, tracer_mass )
333  call thermodyn_r ( rtot, q, tracer_r, qdry )
334  call thermodyn_cp( cptot, q, tracer_cp, qdry )
335  rhot_t_sf(i,j) = sflx_sh(i,j) * rcdz(ks) / ( cptot * gsqrt(ks,i,j,i_xyz) ) &
336  * rhot(ks,i,j) * rtot / pres(ks,i,j) ! = POTT/TEMP
337  enddo
338  enddo
339 
340  if ( i_qv > 0 ) then
341  do j = js, je
342  do i = is, ie
343  work = sflx_qtrc(i,j,i_qv) * rcdz(ks) / gsqrt(ks,i,j,i_xyz)
344  dens_t_sf(i,j) = work
345  rhoq_t_sf(i,j,i_qv) = work
346  enddo
347  enddo
348 
349  do j = js, je
350  do i = is, ie
351  rhot_t_sf(i,j) = rhot_t_sf(i,j) + dens_t_sf(i,j) * rhot(ks,i,j) / dens(ks,i,j)
352  enddo
353  enddo
354 
355  else
356  dens_t_sf(:,:) = 0.0_rp
357  end if
358 
359  endif
360 
361  do j = js, je
362  do i = is, ie
363  dens_t(ks,i,j) = dens_t(ks,i,j) + dens_t_sf(i,j)
364  momz_t(ks,i,j) = momz_t(ks,i,j) + momz_t_sf(i,j)
365  momx_t(ks,i,j) = momx_t(ks,i,j) + momx_t_sf(i,j)
366  momy_t(ks,i,j) = momy_t(ks,i,j) + momy_t_sf(i,j)
367  rhot_t(ks,i,j) = rhot_t(ks,i,j) + rhot_t_sf(i,j)
368  enddo
369  enddo
370 
371  if ( i_qv > 0 ) then
372  do j = js, je
373  do i = is, ie
374  rhoq_t(ks,i,j,i_qv) = rhoq_t(ks,i,j,i_qv) + rhoq_t_sf(i,j,i_qv)
375  enddo
376  enddo
377  end if
378 
379  if ( statistics_checktotal ) then
380  call stat_total( total, dens_t_sf(:,:), 'DENS_t_SF' )
381  call stat_total( total, momz_t_sf(:,:), 'MOMZ_t_SF' )
382  call stat_total( total, momx_t_sf(:,:), 'MOMX_t_SF' )
383  call stat_total( total, momy_t_sf(:,:), 'MOMY_t_SF' )
384  call stat_total( total, rhot_t_sf(:,:), 'RHOT_t_SF' )
385 
386  if ( i_qv > 0 ) then
387  call stat_total( total, rhoq_t_sf(:,:,i_qv), trim(tracer_name(i_qv))//'_t_SF' )
388  end if
389  endif
390 
391  return
real(rp), dimension(:,:,:), allocatable, public dens_tp
logical, public statistics_checktotal
calc&report variable totals to logfile?
integer, public i_xvz
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momz_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sflx_qtrc
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_rhoq_t
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:,:), allocatable, public momy_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_v10
module ATMOSPHERE / Bottom boundary treatment
real(rp), dimension(:,:,:), allocatable, public pres
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhot_t
real(rp), dimension(:,:,:), allocatable, public rhot_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
real(rp), dimension(:), allocatable, public grid_rcdz
reciprocal of center-dz
integer, parameter, public i_lw
module Statistics
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_dn
module Atmosphere / Physics Radiation
module ATMOSPHERIC Surface Variables
integer, parameter, public i_sw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10
real(rp), dimension(:,:,:), allocatable, public temp
module GRIDTRANS
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momx_t
real(rp), dimension(:,:), allocatable, public real_z1
Height of the lowermost grid from surface (cell center) [m].
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_dens_t
module GRID (real space)
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
procedure(sf), pointer, public atmos_phy_sf
integer, public i_xyw
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
module COMMUNICATION
Definition: scale_comm.F90:23
module ATMOSPHERE / Hydrostatic barance
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
module TIME
Definition: scale_time.F90:15
real(dp), public time_dtsec_atmos_phy_sf
time interval of physics(surface flux) [sec]
Definition: scale_time.F90:43
module ATMOSPHERE / Physics Surface fluxes
real(rp), dimension(:,:,:), pointer, public dens_av
logical, public cpl_sw
subroutine, public atmos_bottom_estimate(DENS, PRES, CZ, Zsfc, Z1, SFC_DENS, SFC_PRES)
Calc bottom boundary of atmosphere (just above surface)
real(rp), dimension(:,:,:), allocatable, public momx_tp
module GRID (cartesian)
real(rp), dimension(:,:,:), allocatable, public momz_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momy_t
real(rp), dimension(:,:,:), allocatable, public v
real(rp), dimension(:,:,:,:), allocatable, public gtrans_gsqrt
transformation metrics from Z to Xi, {G}^1/2
real(rp), dimension(:,:,:), allocatable, public w
module ATMOSPHERE / Thermodynamics
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
integer, public i_uyz
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
real(rp), dimension(:,:,:), allocatable, public u
module HISTORY
module TOPOGRAPHY
real(rp), dimension(:), allocatable, public grid_rfdz
reciprocal of face-dz
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
real(rp), dimension(:,:,:), pointer, public rhot_av
integer, public i_xyz
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
module Coupler admin
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Here is the call graph for this function:
Here is the caller graph for this function: