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
  • PARAM_ATMOS_PHY_SF
    nametypedefault valuecomment
    ATMOS_PHY_SF_BETA real(RP) 1.0_RP

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/m2/s SFLX_MU
MVFLX v-momentum flux (merged) kg/m2/s SFLX_MV
MWFLX w-momentum flux (merged) kg/m2/s 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 58 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_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, and scale_process::prc_mpistop().

Referenced by mod_atmos_driver::atmos_driver_setup().

58  use scale_process, only: &
60  use scale_atmos_phy_sf, only: &
62  use mod_atmos_admin, only: &
65  use mod_atmos_phy_sf_vars, only: &
66  sfc_z0m => atmos_phy_sf_sfc_z0m, &
67  sfc_z0h => atmos_phy_sf_sfc_z0h, &
68  sfc_z0e => atmos_phy_sf_sfc_z0e, &
69  sflx_mw => atmos_phy_sf_sflx_mw, &
70  sflx_mu => atmos_phy_sf_sflx_mu, &
71  sflx_mv => atmos_phy_sf_sflx_mv, &
72  sflx_sh => atmos_phy_sf_sflx_sh, &
73  sflx_lh => atmos_phy_sf_sflx_lh, &
74  sflx_qtrc => atmos_phy_sf_sflx_qtrc
75  use mod_cpl_admin, only: &
76  cpl_sw
77  implicit none
78 
79  namelist / param_atmos_phy_sf / &
80  atmos_phy_sf_beta
81 
82  integer :: ierr
83  !---------------------------------------------------------------------------
84 
85  if( io_l ) write(io_fid_log,*)
86  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[ATMOS PHY_SF] / Origin[SCALE-RM]'
87 
88  !--- read namelist
89  rewind(io_fid_conf)
90  read(io_fid_conf,nml=param_atmos_phy_sf,iostat=ierr)
91  if( ierr < 0 ) then !--- missing
92  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
93  elseif( ierr > 0 ) then !--- fatal error
94  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF. Check!'
95  call prc_mpistop
96  endif
97  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_sf)
98 
99  if ( atmos_sw_phy_sf ) then
100 
101  ! setup library component
103 
104  if ( .NOT. cpl_sw ) then
105  if( io_l ) write(io_fid_log,*) '*** Coupler is disabled.'
106  endif
107 
108  else
109 
110  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
111  if( io_l ) write(io_fid_log,*) '*** surface fluxes are set to zero.'
112  sflx_mw(:,:) = 0.0_rp
113  sflx_mu(:,:) = 0.0_rp
114  sflx_mv(:,:) = 0.0_rp
115  sflx_sh(:,:) = 0.0_rp
116  sflx_lh(:,:) = 0.0_rp
117  sflx_qtrc(:,:,:) = 0.0_rp
118 
119  if( io_l ) write(io_fid_log,*) '*** SFC_TEMP, SFC_albedo is set in ATMOS_PHY_SF_vars.'
120 
121  endif
122 
123  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
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
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
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
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 129 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().

129  use mod_atmos_admin, only: &
131  implicit none
132 
133  if ( atmos_sw_phy_sf ) then
134 
135  ! run once (only for the diagnostic value)
136  call prof_rapstart('ATM_SurfaceFlux', 1)
137  call atmos_phy_sf_driver( update_flag = .true. )
138  call prof_rapend ('ATM_SurfaceFlux', 1)
139 
140  end if
141 
142  return
module ATMOS admin
logical, public atmos_sw_phy_sf
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
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 148 of file mod_atmos_phy_sf_driver.f90.

References scale_tracer::aq_name, 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, scale_const::const_cpdry, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_rvap, mod_cpl_admin::cpl_sw, mod_atmos_vars::dens, mod_atmos_vars::dens_tp, scale_grid::grid_rcdz, scale_grid::grid_rfdz, scale_gridtrans::gtrans_gsqrt, scale_tracer::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, mod_atmos_vars::qtrc, scale_grid_real::real_cz, scale_grid_real::real_z1, mod_atmos_vars::rhoq_tp, mod_atmos_vars::rhot, 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, 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().

148  use scale_const, only: &
149  pre00 => const_pre00, &
150  rdry => const_rdry, &
151  rvap => const_rvap, &
152  cpdry => const_cpdry
153  use scale_comm, only: &
154  comm_vars8, &
155  comm_wait
156  use scale_grid, only: &
157  rcdz => grid_rcdz, &
158  rfdz => grid_rfdz
159  use scale_grid_real, only: &
160  real_cz, &
161  real_z1
162  use scale_gridtrans, only: &
163  gsqrt => gtrans_gsqrt, &
164  i_xyz, &
165  i_uyz, &
166  i_xvz, &
167  i_xyw
168  use scale_topography, only: &
169  topo_zsfc
170  use scale_time, only: &
171  dt_sf => time_dtsec_atmos_phy_sf
172  use scale_rm_statistics, only: &
174  stat_total
175  use scale_history, only: &
176  hist_in
177  use scale_atmos_bottom, only: &
178  bottom_estimate => atmos_bottom_estimate
179  use scale_atmos_hydrostatic, only: &
180  barometric_law_mslp => atmos_hydrostatic_barometric_law_mslp
181  use scale_atmos_thermodyn, only: &
182  thermodyn_qd => atmos_thermodyn_qd, &
183  thermodyn_cp => atmos_thermodyn_cp
184  use scale_atmos_phy_sf, only: &
186  use mod_atmos_vars, only: &
187  dens, &
188  rhot, &
189  qtrc, &
190  temp, &
191  pres, &
192  w, &
193  u, &
194  v, &
195  dens_t => dens_tp, &
196  momz_t => momz_tp, &
197  momx_t => momx_tp, &
198  momy_t => momy_tp, &
199  rhot_t => rhot_tp, &
200  rhoq_t => rhoq_tp
201  use mod_atmos_phy_rd_vars, only: &
202  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
203  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn
204  use mod_atmos_phy_sf_vars, only: &
205  dens_t_sf => atmos_phy_sf_dens_t, &
206  momz_t_sf => atmos_phy_sf_momz_t, &
207  momx_t_sf => atmos_phy_sf_momx_t, &
208  momy_t_sf => atmos_phy_sf_momy_t, &
209  rhot_t_sf => atmos_phy_sf_rhot_t, &
210  rhoq_t_sf => atmos_phy_sf_rhoq_t, &
211  sfc_dens => atmos_phy_sf_sfc_dens, &
212  sfc_pres => atmos_phy_sf_sfc_pres, &
213  sfc_temp => atmos_phy_sf_sfc_temp, &
214  sfc_albedo => atmos_phy_sf_sfc_albedo, &
215  sfc_z0m => atmos_phy_sf_sfc_z0m, &
216  sfc_z0h => atmos_phy_sf_sfc_z0h, &
217  sfc_z0e => atmos_phy_sf_sfc_z0e, &
218  sflx_mw => atmos_phy_sf_sflx_mw, &
219  sflx_mu => atmos_phy_sf_sflx_mu, &
220  sflx_mv => atmos_phy_sf_sflx_mv, &
221  sflx_sh => atmos_phy_sf_sflx_sh, &
222  sflx_lh => atmos_phy_sf_sflx_lh, &
223  sflx_gh => atmos_phy_sf_sflx_gh, &
224  sflx_qtrc => atmos_phy_sf_sflx_qtrc, &
225  u10 => atmos_phy_sf_u10, &
226  v10 => atmos_phy_sf_v10, &
227  t2 => atmos_phy_sf_t2, &
228  q2 => atmos_phy_sf_q2
229  use mod_cpl_admin, only: &
230  cpl_sw
231  implicit none
232 
233  logical, intent(in) :: update_flag
234 
235  real(RP) :: uabs10(ia,ja) ! 10m absolute wind [m/s]
236  real(RP) :: mslp (ia,ja) ! mean sea-level pressure [Pa]
237  real(RP) :: beta (ia,ja)
238  real(RP) :: total ! dummy
239 
240  real(RP) :: q(qa)
241  real(RP) :: qdry
242  real(RP) :: rtot
243  real(RP) :: cptot
244 
245  real(RP) :: work
246 
247  integer :: i, j, iq
248  !---------------------------------------------------------------------------
249 
250  if ( update_flag ) then
251 
252  ! update surface density, surface pressure
253  call bottom_estimate( dens(:,:,:), & ! [IN]
254  pres(:,:,:), & ! [IN]
255  real_cz(:,:,:), & ! [IN]
256  topo_zsfc(:,:), & ! [IN]
257  real_z1(:,:), & ! [IN]
258  sfc_dens(:,:), & ! [OUT]
259  sfc_pres(:,:) ) ! [OUT]
260 
261  if ( .NOT. cpl_sw ) then
262  beta(:,:) = atmos_phy_sf_beta
263 
264  call atmos_phy_sf( temp(ks,:,:), & ! [IN]
265  pres(ks,:,:), & ! [IN]
266  w(ks,:,:), & ! [IN]
267  u(ks,:,:), & ! [IN]
268  v(ks,:,:), & ! [IN]
269  dens(ks,:,:), & ! [IN]
270  qtrc(ks,:,:,:), & ! [IN]
271  real_z1(:,:), & ! [IN]
272  dt_sf, & ! [IN]
273  sfc_dens(:,:), & ! [IN]
274  sfc_pres(:,:), & ! [IN]
275  sflx_lw_dn(:,:), & ! [IN]
276  sflx_sw_dn(:,:), & ! [IN]
277  sfc_temp(:,:), & ! [IN]
278  sfc_albedo(:,:,:), & ! [IN]
279  beta(:,:), & ! [IN]
280  sfc_z0m(:,:), & ! [INOUT]
281  sfc_z0h(:,:), & ! [INOUT]
282  sfc_z0e(:,:), & ! [INOUT]
283  sflx_mw(:,:), & ! [OUT]
284  sflx_mu(:,:), & ! [OUT]
285  sflx_mv(:,:), & ! [OUT]
286  sflx_sh(:,:), & ! [OUT]
287  sflx_lh(:,:), & ! [OUT]
288  sflx_qtrc(:,:,:), & ! [OUT]
289  u10(:,:), & ! [OUT]
290  v10(:,:), & ! [OUT]
291  t2(:,:), & ! [OUT]
292  q2(:,:) ) ! [OUT]
293  endif
294 
295 !OCL XFILL
296  do j = js, je
297  do i = is, ie
298  uabs10(i,j) = sqrt( u10(i,j)**2 + v10(i,j)**2 )
299  end do
300  end do
301 
302  call barometric_law_mslp( mslp(:,:), sfc_pres(:,:), t2(:,:), topo_zsfc(:,:) )
303 
304  call hist_in( sfc_dens(:,:), 'SFC_DENS', 'surface atmospheric density', 'kg/m3' )
305  call hist_in( sfc_pres(:,:), 'SFC_PRES', 'surface atmospheric pressure', 'Pa' )
306  call hist_in( sfc_temp(:,:), 'SFC_TEMP', 'surface skin temperature (merged)', 'K' )
307  call hist_in( sfc_albedo(:,:,i_lw), 'SFC_ALB_LW', 'surface albedo (longwave,merged)', '0-1' , nohalo=.true. )
308  call hist_in( sfc_albedo(:,:,i_sw), 'SFC_ALB_SW', 'surface albedo (shortwave,merged)', '0-1' , nohalo=.true. )
309  call hist_in( sfc_z0m(:,:), 'SFC_Z0M', 'roughness length (momentum)', 'm' , nohalo=.true. )
310  call hist_in( sfc_z0h(:,:), 'SFC_Z0H', 'roughness length (heat)', 'm' , nohalo=.true. )
311  call hist_in( sfc_z0e(:,:), 'SFC_Z0E', 'roughness length (vapor)', 'm' , nohalo=.true. )
312  call hist_in( sflx_mw(:,:), 'MWFLX', 'w-momentum flux (merged)', 'kg/m2/s' )
313  call hist_in( sflx_mu(:,:), 'MUFLX', 'u-momentum flux (merged)', 'kg/m2/s' )
314  call hist_in( sflx_mv(:,:), 'MVFLX', 'v-momentum flux (merged)', 'kg/m2/s' )
315  call hist_in( sflx_sh(:,:), 'SHFLX', 'sensible heat flux (merged)', 'W/m2' , nohalo=.true. )
316  call hist_in( sflx_lh(:,:), 'LHFLX', 'latent heat flux (merged)', 'W/m2' , nohalo=.true. )
317  call hist_in( sflx_gh(:,:), 'GHFLX', 'ground heat flux (merged)', 'W/m2' , nohalo=.true. )
318  call hist_in( uabs10(:,:), 'Uabs10', '10m absolute wind', 'm/s' , nohalo=.true. )
319  call hist_in( u10(:,:), 'U10', '10m x-wind', 'm/s' , nohalo=.true. )
320  call hist_in( v10(:,:), 'V10', '10m y-wind', 'm/s' , nohalo=.true. )
321  call hist_in( t2(:,:), 'T2 ', '2m air temperature', 'K' , nohalo=.true. )
322  call hist_in( q2(:,:), 'Q2 ', '2m specific humidity', 'kg/kg' , nohalo=.true. )
323  call hist_in( mslp(:,:), 'MSLP', 'mean sea-level pressure', 'Pa' )
324 
325  call comm_vars8( sflx_mu(:,:), 1 )
326  call comm_vars8( sflx_mv(:,:), 2 )
327  call comm_wait ( sflx_mu(:,:), 1 )
328  call comm_wait ( sflx_mv(:,:), 2 )
329 
330 !OCL XFILL
331  do j = js, je
332  do i = is, ie
333  momz_t_sf(i,j) = sflx_mw(i,j) * rfdz(ks) / gsqrt(ks,i,j,i_xyw)
334  enddo
335  enddo
336 
337 !OCL XFILL
338  do j = js, je
339  do i = is, ie
340  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)
341  enddo
342  enddo
343 
344 !OCL XFILL
345  do j = js, je
346  do i = is, ie
347  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)
348  enddo
349  enddo
350 
351  do j = js, je
352  do i = is, ie
353  q(:) = qtrc(ks,i,j,:)
354  call thermodyn_qd( qdry, q )
355  call thermodyn_cp( cptot, q, qdry )
356  rtot = rdry * qdry + rvap * q(i_qv)
357  rhot_t_sf(i,j) = ( sflx_sh(i,j) * rcdz(ks) / ( cptot * gsqrt(ks,i,j,i_xyz) ) ) &
358  * rhot(ks,i,j) * rtot / pres(ks,i,j) ! = POTT/TEMP
359  enddo
360  enddo
361 
362  dens_t_sf(:,:) = 0.0_rp
363  do iq = i_qv, i_qv
364  do j = js, je
365  do i = is, ie
366  work = sflx_qtrc(i,j,iq) * rcdz(ks) / gsqrt(ks,i,j,i_xyz)
367  dens_t_sf(i,j) = dens_t_sf(i,j) + work
368  rhoq_t_sf(i,j,iq) = work
369  enddo
370  enddo
371  enddo
372 
373  do j = js, je
374  do i = is, ie
375  rhot_t_sf(i,j) = rhot_t_sf(i,j) + dens_t_sf(i,j) * rhot(ks,i,j) / dens(ks,i,j)
376  enddo
377  enddo
378 
379  endif
380 
381  do j = js, je
382  do i = is, ie
383  dens_t(ks,i,j) = dens_t(ks,i,j) + dens_t_sf(i,j)
384  momz_t(ks,i,j) = momz_t(ks,i,j) + momz_t_sf(i,j)
385  momx_t(ks,i,j) = momx_t(ks,i,j) + momx_t_sf(i,j)
386  momy_t(ks,i,j) = momy_t(ks,i,j) + momy_t_sf(i,j)
387  rhot_t(ks,i,j) = rhot_t(ks,i,j) + rhot_t_sf(i,j)
388  enddo
389  enddo
390 
391  do iq = i_qv, i_qv
392  do j = js, je
393  do i = is, ie
394  rhoq_t(ks,i,j,iq) = rhoq_t(ks,i,j,iq) + rhoq_t_sf(i,j,iq)
395  enddo
396  enddo
397  enddo
398 
399  if ( statistics_checktotal ) then
400  call stat_total( total, dens_t_sf(:,:), 'DENS_t_SF' )
401  call stat_total( total, momz_t_sf(:,:), 'MOMZ_t_SF' )
402  call stat_total( total, momx_t_sf(:,:), 'MOMX_t_SF' )
403  call stat_total( total, momy_t_sf(:,:), 'MOMY_t_SF' )
404  call stat_total( total, rhot_t_sf(:,:), 'RHOT_t_SF' )
405 
406  do iq = i_qv, i_qv
407  call stat_total( total, rhoq_t_sf(:,:,iq), trim(aq_name(iq))//'_t_SF' )
408  enddo
409  endif
410 
411  return
real(rp), dimension(:,:,:), allocatable, public dens_tp
integer, public is
start point of inner domain: x, local
logical, public statistics_checktotal
calc&report variable totals to logfile?
integer, public i_xvz
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
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(:,:), allocatable, public atmos_phy_sf_rhot_t
real(rp), dimension(:,:,:), allocatable, public rhot_tp
integer, public qa
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
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
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
integer, public ia
of x whole cells (local, with HALO)
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
integer, public i_qv
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
procedure(sf), pointer, public atmos_phy_sf
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:93
integer, public i_xyw
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
character(len=h_short), dimension(:), allocatable, public aq_name
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
module COMMUNICATION
Definition: scale_comm.F90:23
module ATMOSPHERE / Hydrostatic barance
integer, public js
start point of inner domain: y, local
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
logical, public cpl_sw
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
subroutine, public atmos_bottom_estimate(DENS, PRES, CZ, Zsfc, Z1, SFC_DENS, SFC_PRES)
Calc bottom boundary of atmosphere (just above surface)
integer, public ks
start point of inner domain: z, local
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
integer, public ie
end point of inner domain: x, local
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
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
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function: