SCALE-RM
mod_atmos_phy_sf_driver.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
12 !-------------------------------------------------------------------------------
13 #include "inc_openmp.h"
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
23  use scale_tracer
24 
25  use scale_const, only: &
26  i_sw => const_i_sw, &
27  i_lw => const_i_lw
28  !-----------------------------------------------------------------------------
29  implicit none
30  private
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public procedure
34  !
37  public :: atmos_phy_sf_driver
38 
39  !-----------------------------------------------------------------------------
40  !
41  !++ Public parameters & variables
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private procedure
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private parameters & variables
50  !
51  !-----------------------------------------------------------------------------
52 contains
53  !-----------------------------------------------------------------------------
55  subroutine atmos_phy_sf_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  integer :: ierr
78  !---------------------------------------------------------------------------
79 
80  if( io_l ) write(io_fid_log,*)
81  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[ATMOS PHY_SF] / Origin[SCALE-RM]'
82 
83  if ( atmos_sw_phy_sf ) then
84 
85  ! setup library component
87 
88  if ( .NOT. cpl_sw ) then
89  if( io_l ) write(io_fid_log,*) '*** Coupler is disabled.'
90  endif
91 
92  else
93 
94  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
95  if( io_l ) write(io_fid_log,*) '*** surface fluxes are set to zero.'
96  sflx_mw(:,:) = 0.0_rp
97  sflx_mu(:,:) = 0.0_rp
98  sflx_mv(:,:) = 0.0_rp
99  sflx_sh(:,:) = 0.0_rp
100  sflx_lh(:,:) = 0.0_rp
101  sflx_qtrc(:,:,:) = 0.0_rp
102 
103  if( io_l ) write(io_fid_log,*) '*** SFC_TEMP, SFC_albedo is set in ATMOS_PHY_SF_vars.'
104 
105  endif
106 
107  return
108  end subroutine atmos_phy_sf_driver_setup
109 
110  !-----------------------------------------------------------------------------
112  subroutine atmos_phy_sf_driver_resume
113  use mod_atmos_admin, only: &
115  implicit none
116 
117  if ( atmos_sw_phy_sf ) then
118 
119  ! run once (only for the diagnostic value)
120  call prof_rapstart('ATM_SurfaceFlux', 1)
121  call atmos_phy_sf_driver( update_flag = .true. )
122  call prof_rapend ('ATM_SurfaceFlux', 1)
123 
124  end if
125 
126  return
127  end subroutine atmos_phy_sf_driver_resume
128 
129  !-----------------------------------------------------------------------------
131  subroutine atmos_phy_sf_driver( update_flag )
132  use scale_const, only: &
133  pre00 => const_pre00, &
134  rdry => const_rdry, &
135  rvap => const_rvap, &
136  cpdry => const_cpdry
137  use scale_comm, only: &
138  comm_vars8, &
139  comm_wait
140  use scale_grid, only: &
141  rcdz => grid_rcdz, &
142  rfdz => grid_rfdz
143  use scale_grid_real, only: &
144  real_cz, &
145  real_z1
146  use scale_gridtrans, only: &
147  gsqrt => gtrans_gsqrt, &
148  i_xyz, &
149  i_uyz, &
150  i_xvz, &
151  i_xyw
152  use scale_topography, only: &
153  topo_zsfc
154  use scale_time, only: &
155  dt_sf => time_dtsec_atmos_phy_sf
156  use scale_rm_statistics, only: &
158  stat_total
159  use scale_history, only: &
160  hist_in
161  use scale_atmos_bottom, only: &
162  bottom_estimate => atmos_bottom_estimate
163  use scale_atmos_hydrostatic, only: &
164  barometric_law_mslp => atmos_hydrostatic_barometric_law_mslp
165  use scale_atmos_thermodyn, only: &
166  thermodyn_qd => atmos_thermodyn_qd, &
167  thermodyn_cp => atmos_thermodyn_cp
168  use scale_atmos_phy_sf, only: &
170  use mod_atmos_vars, only: &
171  dens, &
172  rhot, &
173  qtrc, &
174  temp, &
175  pres, &
176  w, &
177  u, &
178  v, &
179  dens_t => dens_tp, &
180  momz_t => momz_tp, &
181  momx_t => momx_tp, &
182  momy_t => momy_tp, &
183  rhot_t => rhot_tp, &
184  rhoq_t => rhoq_tp
185  use mod_atmos_phy_rd_vars, only: &
186  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
187  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn
188  use mod_atmos_phy_sf_vars, only: &
189  dens_t_sf => atmos_phy_sf_dens_t, &
190  momz_t_sf => atmos_phy_sf_momz_t, &
191  momx_t_sf => atmos_phy_sf_momx_t, &
192  momy_t_sf => atmos_phy_sf_momy_t, &
193  rhot_t_sf => atmos_phy_sf_rhot_t, &
194  rhoq_t_sf => atmos_phy_sf_rhoq_t, &
195  sfc_dens => atmos_phy_sf_sfc_dens, &
196  sfc_pres => atmos_phy_sf_sfc_pres, &
197  sfc_temp => atmos_phy_sf_sfc_temp, &
198  sfc_albedo => atmos_phy_sf_sfc_albedo, &
199  sfc_z0m => atmos_phy_sf_sfc_z0m, &
200  sfc_z0h => atmos_phy_sf_sfc_z0h, &
201  sfc_z0e => atmos_phy_sf_sfc_z0e, &
202  sflx_mw => atmos_phy_sf_sflx_mw, &
203  sflx_mu => atmos_phy_sf_sflx_mu, &
204  sflx_mv => atmos_phy_sf_sflx_mv, &
205  sflx_sh => atmos_phy_sf_sflx_sh, &
206  sflx_lh => atmos_phy_sf_sflx_lh, &
207  sflx_gh => atmos_phy_sf_sflx_gh, &
208  sflx_qtrc => atmos_phy_sf_sflx_qtrc, &
209  u10 => atmos_phy_sf_u10, &
210  v10 => atmos_phy_sf_v10, &
211  t2 => atmos_phy_sf_t2, &
212  q2 => atmos_phy_sf_q2
213  use mod_cpl_admin, only: &
214  cpl_sw
215  implicit none
216 
217  logical, intent(in) :: update_flag
218 
219  real(RP) :: Uabs10(ia,ja) ! 10m absolute wind [m/s]
220  real(RP) :: MSLP (ia,ja) ! mean sea-level pressure [Pa]
221  real(RP) :: total ! dummy
222 
223  real(RP) :: q(qa)
224  real(RP) :: qdry
225  real(RP) :: Rtot
226  real(RP) :: CPtot
227 
228  real(RP) :: work
229 
230  integer :: i, j, iq
231  !---------------------------------------------------------------------------
232 
233  if ( update_flag ) then
234 
235  ! update surface density, surface pressure
236  call bottom_estimate( dens(:,:,:), & ! [IN]
237  pres(:,:,:), & ! [IN]
238  real_cz(:,:,:), & ! [IN]
239  topo_zsfc(:,:), & ! [IN]
240  real_z1(:,:), & ! [IN]
241  sfc_dens(:,:), & ! [OUT]
242  sfc_pres(:,:) ) ! [OUT]
243 
244  if ( .NOT. cpl_sw ) then
245  call atmos_phy_sf( temp(ks,:,:), & ! [IN]
246  pres(ks,:,:), & ! [IN]
247  w(ks,:,:), & ! [IN]
248  u(ks,:,:), & ! [IN]
249  v(ks,:,:), & ! [IN]
250  dens(ks,:,:), & ! [IN]
251  qtrc(ks,:,:,:), & ! [IN]
252  real_z1(:,:), & ! [IN]
253  dt_sf, & ! [IN]
254  sfc_dens(:,:), & ! [IN]
255  sfc_pres(:,:), & ! [IN]
256  sflx_lw_dn(:,:), & ! [IN]
257  sflx_sw_dn(:,:), & ! [IN]
258  sfc_temp(:,:), & ! [IN]
259  sfc_albedo(:,:,:), & ! [IN]
260  sfc_z0m(:,:), & ! [INOUT]
261  sfc_z0h(:,:), & ! [INOUT]
262  sfc_z0e(:,:), & ! [INOUT]
263  sflx_mw(:,:), & ! [OUT]
264  sflx_mu(:,:), & ! [OUT]
265  sflx_mv(:,:), & ! [OUT]
266  sflx_sh(:,:), & ! [OUT]
267  sflx_lh(:,:), & ! [OUT]
268  sflx_qtrc(:,:,:), & ! [OUT]
269  u10(:,:), & ! [OUT]
270  v10(:,:), & ! [OUT]
271  t2(:,:), & ! [OUT]
272  q2(:,:) ) ! [OUT]
273  endif
274 
275 !OCL XFILL
276  do j = js, je
277  do i = is, ie
278  uabs10(i,j) = sqrt( u10(i,j)**2 + v10(i,j)**2 )
279  end do
280  end do
281 
282  call barometric_law_mslp( mslp(:,:), sfc_pres(:,:), t2(:,:), topo_zsfc(:,:) )
283 
284  call hist_in( sfc_dens(:,:), 'SFC_DENS', 'surface atmospheric density', 'kg/m3' )
285  call hist_in( sfc_pres(:,:), 'SFC_PRES', 'surface atmospheric pressure', 'Pa' )
286  call hist_in( sfc_temp(:,:), 'SFC_TEMP', 'surface skin temperature (merged)', 'K' )
287  call hist_in( sfc_albedo(:,:,i_lw), 'SFC_ALB_LW', 'surface albedo (longwave,merged)', '0-1' , nohalo=.true. )
288  call hist_in( sfc_albedo(:,:,i_sw), 'SFC_ALB_SW', 'surface albedo (shortwave,merged)', '0-1' , nohalo=.true. )
289  call hist_in( sfc_z0m(:,:), 'SFC_Z0M', 'roughness length (momentum)', 'm' , nohalo=.true. )
290  call hist_in( sfc_z0h(:,:), 'SFC_Z0H', 'roughness length (heat)', 'm' , nohalo=.true. )
291  call hist_in( sfc_z0e(:,:), 'SFC_Z0E', 'roughness length (vapor)', 'm' , nohalo=.true. )
292  call hist_in( sflx_mw(:,:), 'MWFLX', 'w-momentum flux (merged)', 'kg/m2/s' )
293  call hist_in( sflx_mu(:,:), 'MUFLX', 'u-momentum flux (merged)', 'kg/m2/s' )
294  call hist_in( sflx_mv(:,:), 'MVFLX', 'v-momentum flux (merged)', 'kg/m2/s' )
295  call hist_in( sflx_sh(:,:), 'SHFLX', 'sensible heat flux (merged)', 'W/m2' , nohalo=.true. )
296  call hist_in( sflx_lh(:,:), 'LHFLX', 'latent heat flux (merged)', 'W/m2' , nohalo=.true. )
297  call hist_in( sflx_gh(:,:), 'GHFLX', 'ground heat flux (merged)', 'W/m2' , nohalo=.true. )
298  call hist_in( uabs10(:,:), 'Uabs10', '10m absolute wind', 'm/s' , nohalo=.true. )
299  call hist_in( u10(:,:), 'U10', '10m x-wind', 'm/s' , nohalo=.true. )
300  call hist_in( v10(:,:), 'V10', '10m y-wind', 'm/s' , nohalo=.true. )
301  call hist_in( t2(:,:), 'T2 ', '2m air temperature', 'K' , nohalo=.true. )
302  call hist_in( q2(:,:), 'Q2 ', '2m specific humidity', 'kg/kg' , nohalo=.true. )
303  call hist_in( mslp(:,:), 'MSLP', 'mean sea-level pressure', 'Pa' )
304 
305  call comm_vars8( sflx_mu(:,:), 1 )
306  call comm_vars8( sflx_mv(:,:), 2 )
307  call comm_wait ( sflx_mu(:,:), 1 )
308  call comm_wait ( sflx_mv(:,:), 2 )
309 
310 !OCL XFILL
311  do j = js, je
312  do i = is, ie
313  momz_t_sf(i,j) = sflx_mw(i,j) * rfdz(ks) / gsqrt(ks,i,j,i_xyw)
314  enddo
315  enddo
316 
317 !OCL XFILL
318  do j = js, je
319  do i = is, ie
320  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)
321  enddo
322  enddo
323 
324 !OCL XFILL
325  do j = js, je
326  do i = is, ie
327  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)
328  enddo
329  enddo
330 
331  do j = js, je
332  do i = is, ie
333  q(:) = qtrc(ks,i,j,:)
334  call thermodyn_qd( qdry, q )
335  call thermodyn_cp( cptot, q, qdry )
336  rtot = rdry * qdry + rvap * q(i_qv)
337  rhot_t_sf(i,j) = ( sflx_sh(i,j) * rcdz(ks) / ( cptot * gsqrt(ks,i,j,i_xyz) ) ) &
338  * rhot(ks,i,j) * rtot / pres(ks,i,j) ! = POTT/TEMP
339  enddo
340  enddo
341 
342  dens_t_sf(:,:) = 0.0_rp
343  do iq = i_qv, i_qv
344  do j = js, je
345  do i = is, ie
346  work = sflx_qtrc(i,j,iq) * rcdz(ks) / gsqrt(ks,i,j,i_xyz)
347  dens_t_sf(i,j) = dens_t_sf(i,j) + work
348  rhoq_t_sf(i,j,iq) = work
349  enddo
350  enddo
351  enddo
352 
353  do j = js, je
354  do i = is, ie
355  rhot_t_sf(i,j) = rhot_t_sf(i,j) + dens_t_sf(i,j) * rhot(ks,i,j) / dens(ks,i,j)
356  enddo
357  enddo
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  do iq = i_qv, i_qv
372  do j = js, je
373  do i = is, ie
374  rhoq_t(ks,i,j,iq) = rhoq_t(ks,i,j,iq) + rhoq_t_sf(i,j,iq)
375  enddo
376  enddo
377  enddo
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  do iq = i_qv, i_qv
387  call stat_total( total, rhoq_t_sf(:,:,iq), trim(aq_name(iq))//'_t_SF' )
388  enddo
389  endif
390 
391  return
392  end subroutine atmos_phy_sf_driver
393 
394 end module mod_atmos_phy_sf_driver
module ATMOS admin
real(rp), dimension(:,:,:), allocatable, public dens_tp
integer, public is
start point of inner domain: x, local
subroutine, public atmos_phy_sf_setup(SF_TYPE)
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
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
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_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
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module ATMOSPHERE / Bottom boundary treatment
real(rp), dimension(:,:,:), allocatable, public pres
subroutine, public atmos_phy_sf_driver_resume
Resume.
module ATMOSPHERIC Variables
module ATMOSPHERE / Physics Surface fluxes
module STDIO
Definition: scale_stdio.F90:12
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
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 grid index
character(len=h_short), public atmos_phy_sf_type
module ATMOSPHERIC Surface Variables
module TRACER
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
subroutine, public atmos_phy_sf_driver_setup
Setup.
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
module PROCESS
real(dp), public time_dtsec_atmos_phy_sf
time interval of physics(surface flux) [sec]
Definition: scale_time.F90:43
logical, public atmos_sw_phy_sf
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)
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
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
module profiler
Definition: scale_prof.F90:10
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
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
module PRECISION
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
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_gh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
subroutine, public atmos_phy_sf_driver(update_flag)
Driver.
integer, public ja
of y whole cells (local, with HALO)