SCALE-RM
mod_atmos_phy_sf_driver.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_tracer
23  !-----------------------------------------------------------------------------
24  implicit none
25  private
26  !-----------------------------------------------------------------------------
27  !
28  !++ Public procedure
29  !
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private procedure
40  !
41  !-----------------------------------------------------------------------------
42  !
43  !++ Private parameters & variables
44  !
45  integer :: hist_uabs10, hist_u10m, hist_v10m
46  integer :: hist_t2, hist_q2, hist_rh2
47  integer :: hist_mslp
48 
49  !-----------------------------------------------------------------------------
50 contains
51  !-----------------------------------------------------------------------------
53  subroutine atmos_phy_sf_driver_setup
54  use scale_prc, only: &
55  prc_abort
56  use scale_atmos_phy_sf_bulk, only: &
58  use scale_atmos_phy_sf_const, only: &
60  use scale_file_history, 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_shex => atmos_phy_sf_sflx_shex, &
75  sflx_qvex => atmos_phy_sf_sflx_qvex, &
76  sflx_qtrc => atmos_phy_sf_sflx_qtrc, &
77  sflx_engi => atmos_phy_sf_sflx_engi, &
78  ustar => atmos_phy_sf_ustar, &
79  tstar => atmos_phy_sf_tstar, &
80  qstar => atmos_phy_sf_qstar, &
81  wstar => atmos_phy_sf_wstar, &
82  rlmo => atmos_phy_sf_rlmo
83  use mod_cpl_admin, only: &
84  cpl_sw
85  implicit none
86  !---------------------------------------------------------------------------
87 
88  log_newline
89  log_info("ATMOS_PHY_SF_driver_setup",*) 'Setup'
90 
91  if ( atmos_sw_phy_sf ) then
92 
93  if ( cpl_sw ) then
94  log_info("ATMOS_PHY_SF_driver_setup",*) 'Coupler is enabled.'
95  else
96  ! setup library component
97  select case( atmos_phy_sf_type )
98  case ( 'BULK' )
100  case ( 'CONST' )
102  case default
103  log_error("ATMOS_PHY_SF_driver_setup",*) 'invalid Surface flux type(', trim(atmos_phy_sf_type), '). CHECK!'
104  call prc_abort
105  end select
106  endif
107 
108  else
109 
110  log_info("ATMOS_PHY_SF_driver_setup",*) 'this component is never called.'
111  log_info("ATMOS_PHY_SF_driver_setup",*) 'surface fluxes are set to zero.'
112  !$acc kernels
113  sflx_mw(:,:) = 0.0_rp
114  sflx_mu(:,:) = 0.0_rp
115  sflx_mv(:,:) = 0.0_rp
116  sflx_sh(:,:) = 0.0_rp
117  sflx_lh(:,:) = 0.0_rp
118  sflx_shex(:,:) = 0.0_rp
119  sflx_qvex(:,:) = 0.0_rp
120  ustar(:,:) = 0.0_rp
121  tstar(:,:) = 0.0_rp
122  qstar(:,:) = 0.0_rp
123  wstar(:,:) = 0.0_rp
124  rlmo(:,:) = 0.0_rp
125  !$acc end kernels
126  log_info("ATMOS_PHY_SF_driver_setup",*) 'SFC_TEMP, SFC_albedo is set in ATMOS_PHY_SF_vars.'
127 
128  endif
129 
130  !$acc kernels
131  sflx_qtrc(:,:,:) = 0.0_rp
132  sflx_engi(:,:) = 0.0_rp
133  !$acc end kernels
134 
135  call file_history_reg( 'Uabs10', '10m absolute wind', 'm/s' , hist_uabs10, ndims=2, fill_halo=.true. )
136  call file_history_reg( 'U10m', '10m eastward wind', 'm/s' , hist_u10m, ndims=2, fill_halo=.true. )
137  call file_history_reg( 'V10m', '10m northward wind', 'm/s' , hist_v10m, ndims=2, fill_halo=.true. )
138  call file_history_reg( 'T2', '2m air temperature', 'K' , hist_t2, ndims=2, fill_halo=.true. )
139  call file_history_reg( 'Q2', '2m specific humidity', 'kg/kg', hist_q2, ndims=2, fill_halo=.true. )
140  call file_history_reg( 'RH2', '2m relative humidity', '%', hist_rh2, ndims=2, fill_halo=.true. )
141  call file_history_reg( 'MSLP', 'mean sea-level pressure', 'Pa' , hist_mslp, ndims=2, fill_halo=.true., standard_name='air_pressure_at_mean_sea_level' )
142 
143  return
144  end subroutine atmos_phy_sf_driver_setup
145 
146  !-----------------------------------------------------------------------------
148  subroutine atmos_phy_sf_driver_calc_tendency( update_flag )
149  use scale_const, only: &
150  undef => const_undef, &
151  pre00 => const_pre00, &
152  rdry => const_rdry, &
153  rvap => const_rvap, &
154  cpdry => const_cpdry, &
155  cpvap => const_cpvap, &
156  epstvap => const_epstvap
157  use scale_atmos_grid_cartesc_real, only: &
163  use scale_topography, only: &
164  tansl_x => topography_tansl_x, &
165  tansl_y => topography_tansl_y
166  use scale_time, only: &
167  dt_sf => time_dtsec_atmos_phy_sf
168  use scale_statistics, only: &
170  statistics_total
171  use scale_atmos_bottom, only: &
172  bottom_estimate => atmos_bottom_estimate
173  use scale_atmos_hydrometeor, only: &
175  lhv, &
176  i_qv
177  use scale_atmos_phy_sf_bulk, only: &
179  use scale_atmos_phy_sf_const, only: &
181  use scale_bulkflux, only: &
182  bulkflux_diagnose_scales
183  use mod_atmos_admin, only: &
185  use mod_atmos_vars, only: &
186  dens => dens_av, &
187  rhot => rhot_av, &
188  pott, &
189  temp, &
190  pres, &
191  w, &
192  u, &
193  v, &
194  qv, &
195  cptot, &
196  cvtot, &
197  dens_t => dens_tp, &
198  momz_t => momz_tp, &
199  rhou_t => rhou_tp, &
200  rhov_t => rhov_tp, &
201  rhoh => rhoh_p, &
202  rhoq_t => rhoq_tp
203  use mod_atmos_phy_rd_vars, only: &
204  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
205  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn
206  use mod_atmos_phy_bl_vars, only: &
207  pbl_zi => atmos_phy_bl_zi
208  use mod_atmos_phy_sf_vars, only: &
209  dens_t_sf => atmos_phy_sf_dens_t, &
210  momz_t_sf => atmos_phy_sf_momz_t, &
211  rhou_t_sf => atmos_phy_sf_rhou_t, &
212  rhov_t_sf => atmos_phy_sf_rhov_t, &
213  rhoh_sf => atmos_phy_sf_rhoh, &
214  rhoq_t_sf => atmos_phy_sf_rhoq_t, &
215  sfc_dens => atmos_phy_sf_sfc_dens, &
216  sfc_pres => atmos_phy_sf_sfc_pres, &
217  sfc_temp => atmos_phy_sf_sfc_temp, &
218  sfc_z0m => atmos_phy_sf_sfc_z0m, &
219  sfc_z0h => atmos_phy_sf_sfc_z0h, &
220  sfc_z0e => atmos_phy_sf_sfc_z0e, &
221  sflx_mw => atmos_phy_sf_sflx_mw, &
222  sflx_mu => atmos_phy_sf_sflx_mu, &
223  sflx_mv => atmos_phy_sf_sflx_mv, &
224  sflx_sh => atmos_phy_sf_sflx_sh, &
225  sflx_lh => atmos_phy_sf_sflx_lh, &
226  sflx_shex => atmos_phy_sf_sflx_shex, &
227  sflx_qvex => atmos_phy_sf_sflx_qvex, &
228  sflx_qtrc => atmos_phy_sf_sflx_qtrc, &
229  sflx_engi => atmos_phy_sf_sflx_engi, &
230  ustar => atmos_phy_sf_ustar, &
231  tstar => atmos_phy_sf_tstar, &
232  qstar => atmos_phy_sf_qstar, &
233  wstar => atmos_phy_sf_wstar, &
234  rlmo => atmos_phy_sf_rlmo, &
235  u10 => atmos_phy_sf_u10, &
236  v10 => atmos_phy_sf_v10, &
237  t2 => atmos_phy_sf_t2, &
238  q2 => atmos_phy_sf_q2
239  use mod_cpl_admin, only: &
240  cpl_sw
241  implicit none
242 
243  logical, intent(in) :: update_flag
244 
245  real(rp) :: atm_w (ia,ja)
246  real(rp) :: atm_u (ia,ja)
247  real(rp) :: atm_v (ia,ja)
248  real(rp) :: atm_dens(ia,ja)
249  real(rp) :: atm_temp(ia,ja)
250  real(rp) :: atm_pres(ia,ja)
251  real(rp) :: atm_qv (ia,ja)
252  real(rp) :: sfc_potv(ia,ja)
253  real(rp) :: sflx_sh2(ia,ja)
254  real(rp) :: sflx_qv (ia,ja)
255  real(rp) :: cp_t, cv_t
256  real(rp) :: engi_t
257  real(rp) :: rdz
258  real(rp) :: work
259  real(rp) :: kappa
260 
261  integer :: i, j, iq
262  !---------------------------------------------------------------------------
263 
264  if ( update_flag ) then
265 
266  !$acc data create(ATM_DENS,ATM_TEMP,ATM_PRES,SFLX_QV)
267 
268  ! update surface density, surface pressure
269  call bottom_estimate( ka, ks, ke, ia, is, ie, ja, js, je, &
270  dens(:,:,:), pres(:,:,:), qv(:,:,:), & ! [IN]
271  sfc_temp(:,:), & ! [IN]
272  fz(:,:,:), & ! [IN]
273  sfc_dens(:,:), sfc_pres(:,:) ) ! [OUT]
274 
275  !$omp parallel do
276  !$acc kernels
277  do j = js, je
278  do i = is, ie
279  atm_dens(i,j) = dens(ks,i,j)
280  atm_temp(i,j) = temp(ks,i,j)
281  atm_pres(i,j) = pres(ks,i,j)
282  end do
283  end do
284  !$acc end kernels
285 
286  if ( cpl_sw ) then
287 
288  !$acc data create(SFLX_SH2)
289 
290  !$omp parallel do
291  !$acc kernels
292  do j = js, je
293  do i = is, ie
294  sflx_sh2(i,j) = sflx_sh(i,j) - sflx_shex(i,j)
295  end do
296  end do
297  !$acc end kernels
298 
299  if ( atmos_hydrometeor_dry ) then
300  !$omp parallel do
301  !$acc kernels
302  do j = js, je
303  do i = is, ie
304  sflx_qv(i,j) = 0.0_rp
305  sfc_potv(i,j) = sfc_temp(i,j) * ( pre00 / sfc_pres(i,j) )**( rdry / cpdry )
306  end do
307  end do
308  !$acc end kernels
309  else
310  !$omp parallel do private(kappa)
311  !$acc kernels
312  do j = js, je
313  do i = is, ie
314  sflx_qv(i,j) = sflx_qtrc(i,j,i_qv) - sflx_qvex(i,j)
315  kappa = ( rdry + ( rvap - rdry ) * qv(ks,i,j) ) &
316  / ( cpdry + ( cpvap - cpdry ) * qv(ks,i,j) )
317  sfc_potv(i,j) = sfc_temp(i,j) * ( pre00 / sfc_pres(i,j) )**kappa &
318  * ( 1.0_rp + epstvap * qv(ks,i,j) )
319  end do
320  end do
321  !$acc end kernels
322  end if
323 
324  call bulkflux_diagnose_scales( ia, is, ie, ja, js, je, &
325  sflx_mw(:,:), sflx_mu(:,:), sflx_mv(:,:), & ! [IN]
326  sflx_sh2(:,:), sflx_qv(:,:), & ! [IN]
327  sfc_dens(:,:), sfc_potv(:,:), pbl_zi(:,:), & ! [IN]
328  ustar(:,:), tstar(:,:), qstar(:,:), & ! [OUT]
329  wstar(:,:), rlmo(:,:) ) ! [OUT]
330  !$acc end data
331 
332  else
333 
334  !$acc data create(ATM_U,ATM_V,ATM_W,ATM_QV)
335 
336  !$omp parallel do
337  !$acc kernels
338  do j = js, je
339  do i = is, ie
340  atm_u(i,j) = u(ks,i,j)
341  atm_v(i,j) = v(ks,i,j)
342  atm_w(i,j) = atm_u(i,j) * tansl_x(i,j) + atm_v(i,j) * tansl_y(i,j)
343  atm_qv(i,j) = qv(ks,i,j)
344  enddo
345  enddo
346  !$acc end kernels
347 
348  select case ( atmos_phy_sf_type )
349  case ( 'BULK' )
350 
351  call atmos_phy_sf_bulk_flux( ia, is, ie, ja, js, je, & ! [IN]
352  atm_w(:,:), atm_u(:,:), atm_v(:,:), & ! [IN]
353  atm_temp(:,:), atm_pres(:,:), atm_qv(:,:), & ! [IN]
354  sfc_dens(:,:), sfc_temp(:,:), sfc_pres(:,:), & ! [IN]
355  sfc_z0m(:,:), sfc_z0h(:,:), sfc_z0e(:,:), & ! [IN]
356  pbl_zi(:,:), z1(:,:), & ! [IN]
357  sflx_mw(:,:), sflx_mu(:,:), sflx_mv(:,:), & ! [OUT]
358  sflx_sh(:,:), sflx_lh(:,:), sflx_qv(:,:), & ! [OUT]
359  ustar(:,:), tstar(:,:), qstar(:,:), & ! [OUT]
360  wstar(:,:), & ! [OUT]
361  rlmo(:,:), & ! [OUT]
362  u10(:,:), v10(:,:), t2(:,:), q2(:,:) ) ! [OUT]
363 
364  case ( 'CONST' )
365 
366  !$acc update host(ATM_W,ATM_U,ATM_V,ATM_TEMP,SFC_DENS)
367  call atmos_phy_sf_const_flux( ia, is, ie, ja, js, je, & ! [IN]
368  atm_w(:,:), atm_u(:,:), atm_v(:,:), sfc_temp(:,:), & ! [IN]
369  z1(:,:), sfc_dens(:,:), & ! [IN]
370  sflx_mw(:,:), sflx_mu(:,:), sflx_mv(:,:), & ! [OUT]
371  sflx_sh(:,:), sflx_lh(:,:), sflx_qv(:,:), & ! [OUT]
372  u10(:,:), v10(:,:) ) ! [OUT]
373  ustar(:,:) = undef
374  tstar(:,:) = undef
375  qstar(:,:) = undef
376  rlmo(:,:) = undef
377  t2(:,:) = atm_temp(:,:)
378  q2(:,:) = atm_qv(:,:)
379  !$acc update device(SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_QV,Ustar,Tstar,Qstar,Wstar,RLmo,U10,V10,T2,Q2)
380 
381  end select
382 
383  if ( .NOT. atmos_hydrometeor_dry ) then
384  !$acc kernels
385  sflx_qtrc(:,:,i_qv) = sflx_qv(:,:)
386  sflx_engi(:,:) = sflx_qv(:,:) * ( tracer_cv(i_qv) * sfc_temp(:,:) + lhv )
387  !$acc end kernels
388  endif
389 
390  !$acc end data
391 
392  endif
393 
394  call history_output
395 
396 !OCL XFILL
397  !$omp parallel do &
398  !$omp private(rdz)
399  !$acc kernels
400  do j = js, je
401  do i = is, ie
402  rdz = 1.0_rp / ( fz(ks,i,j) - fz(ks-1,i,j) )
403  momz_t_sf(i,j) = sflx_mw(i,j) / ( cz(ks+1,i,j) - cz(ks,i,j) )
404  rhou_t_sf(i,j) = sflx_mu(i,j) * rdz
405  rhov_t_sf(i,j) = sflx_mv(i,j) * rdz
406  enddo
407  enddo
408  !$acc end kernels
409 
410  !$omp parallel do &
411  !$omp private(work,rdz,CP_t,CV_t,ENGI_t)
412  !$acc kernels
413  do j = js, je
414  do i = is, ie
415  rdz = 1.0_rp / ( fz(ks,i,j) - fz(ks-1,i,j) )
416  dens_t_sf(i,j) = 0.0_rp
417  cp_t = 0.0_rp
418  cv_t = 0.0_rp
419  engi_t = sflx_engi(i,j) * rdz
420  !$acc loop seq
421  do iq = 1, qa
422  work = sflx_qtrc(i,j,iq) * rdz
423 
424  rhoq_t_sf(i,j,iq) = work
425  dens_t_sf(i,j) = dens_t_sf(i,j) + work * tracer_mass(iq)
426  cp_t = cp_t + work * tracer_cp(iq)
427  cv_t = cv_t + work * tracer_cv(iq)
428  engi_t = engi_t - tracer_engi0(iq) * work
429  end do
430  cp_t = ( cp_t - cptot(ks,i,j) * dens_t_sf(i,j) ) / atm_dens(i,j)
431  cv_t = ( cv_t - cvtot(ks,i,j) * dens_t_sf(i,j) ) / atm_dens(i,j)
432 
433  rhoh_sf(i,j) = sflx_sh(i,j) * rdz + engi_t &
434  - ( cp_t + log( atm_pres(i,j) / pre00 ) * ( cvtot(ks,i,j) / cptot(ks,i,j) * cp_t - cv_t ) ) * atm_dens(i,j) * atm_temp(i,j)
435  enddo
436  enddo
437  !$acc end kernels
438 
439  !$acc end data
440 
441  endif
442 
443  !$omp parallel do
444  !$acc kernels
445  do j = js, je
446  do i = is, ie
447  momz_t(ks,i,j) = momz_t(ks,i,j) + momz_t_sf(i,j)
448  rhou_t(ks,i,j) = rhou_t(ks,i,j) + rhou_t_sf(i,j)
449  rhov_t(ks,i,j) = rhov_t(ks,i,j) + rhov_t_sf(i,j)
450  rhoh(ks,i,j) = rhoh(ks,i,j) + rhoh_sf(i,j)
451  dens_t(ks,i,j) = dens_t(ks,i,j) + dens_t_sf(i,j)
452  enddo
453  enddo
454  !$acc end kernels
455 
456  !$omp parallel
457  !$acc kernels
458  do iq = 1, qa
459  !$omp do
460  do j = js, je
461  do i = is, ie
462  rhoq_t(ks,i,j,iq) = rhoq_t(ks,i,j,iq) + rhoq_t_sf(i,j,iq)
463  enddo
464  enddo
465  !$omp end do nowait
466  enddo
467  !$acc end kernels
468  !$omp end parallel
469 
470  if ( statistics_checktotal ) then
471 
472  if ( .NOT. atmos_hydrometeor_dry ) then
473  do iq = 1, qa
474  call statistics_total( ia, is, ie, ja, js, je, &
475  sflx_qtrc(:,:,iq), 'SFLX_'//trim(tracer_name(iq)), &
478  enddo
479  endif
480 
481  call statistics_total( ia, is, ie, ja, js, je, &
482  dens_t_sf(:,:), 'DENS_t_SF', &
485  call statistics_total( ia, is, ie, ja, js, je, &
486  momz_t_sf(:,:), 'MOMZ_t_SF', &
489  call statistics_total( ia, is, ie, ja, js, je, &
490  rhou_t_sf(:,:), 'RHOU_t_SF', &
493  call statistics_total( ia, is, ie, ja, js, je, &
494  rhov_t_sf(:,:), 'RHOV_t_SF', &
497  call statistics_total( ia, is, ie, ja, js, je, &
498  rhoh_sf(:,:), 'RHOH_SF', &
501 
502  do iq = 1, qa
503  call statistics_total( ia, is, ie, ja, js, je, &
504  rhoq_t_sf(:,:,iq), trim(tracer_name(iq))//'_t_SF', &
507  enddo
508 
509  endif
510 
511  return
512  end subroutine atmos_phy_sf_driver_calc_tendency
513 
514  !-----------------------------------------------------------------------------
515  subroutine history_output
516  use scale_const, only: &
517  undef => const_undef, &
518  rvap => const_rvap
519  use scale_file_history, only: &
520  file_history_query, &
521  file_history_in, &
522  file_history_put
523  use scale_atmos_hydrostatic, only: &
524  barometric_law_mslp => atmos_hydrostatic_barometric_law_mslp
525  use scale_atmos_saturation, only: &
526  atmos_saturation_psat_liq
527  use scale_atmos_grid_cartesc_real, only: &
528  real_cz => atmos_grid_cartesc_real_cz
531  use scale_atmos_hydrometeor, only: &
533  i_qv
534  use mod_atmos_vars, only: &
535  temp, &
536  pres, &
537  qv
538  use mod_atmos_phy_sf_vars, only: &
539  sfc_dens => atmos_phy_sf_sfc_dens, &
540  sfc_pres => atmos_phy_sf_sfc_pres, &
541  sfc_temp => atmos_phy_sf_sfc_temp, &
542  sfc_albedo => atmos_phy_sf_sfc_albedo, &
543  sfc_z0m => atmos_phy_sf_sfc_z0m, &
544  sfc_z0h => atmos_phy_sf_sfc_z0h, &
545  sfc_z0e => atmos_phy_sf_sfc_z0e, &
546  sflx_mw => atmos_phy_sf_sflx_mw, &
547  sflx_mu => atmos_phy_sf_sflx_mu, &
548  sflx_mv => atmos_phy_sf_sflx_mv, &
549  sflx_sh => atmos_phy_sf_sflx_sh, &
550  sflx_lh => atmos_phy_sf_sflx_lh, &
551  sflx_gh => atmos_phy_sf_sflx_gh, &
552  sflx_qtrc => atmos_phy_sf_sflx_qtrc, &
553  sflx_engi => atmos_phy_sf_sflx_engi, &
554  ustar => atmos_phy_sf_ustar, &
555  tstar => atmos_phy_sf_tstar, &
556  qstar => atmos_phy_sf_qstar, &
557  wstar => atmos_phy_sf_wstar, &
558  rlmo => atmos_phy_sf_rlmo, &
559  u10 => atmos_phy_sf_u10, &
560  v10 => atmos_phy_sf_v10, &
561  t2 => atmos_phy_sf_t2, &
562  q2 => atmos_phy_sf_q2
563  implicit none
564 
565  real(RP) :: work(IA,JA)
566 
567  logical :: do_put
568  integer :: i, j, iq
569  !---------------------------------------------------------------------------
570 
571  call file_history_in( sfc_dens(:,:), 'SFC_DENS', 'surface atmospheric density', 'kg/m3' )
572  call file_history_in( sfc_pres(:,:), 'SFC_PRES', 'surface atmospheric pressure', 'Pa' )
573  call file_history_in( sfc_temp(:,:), 'SFC_TEMP', 'surface skin temperature', 'K' )
574  call file_history_in( sfc_albedo(:,:,i_r_direct ,i_r_ir ), 'SFC_ALB_IR_dir' , 'surface albedo (IR; direct', '1' , fill_halo=.true. )
575  call file_history_in( sfc_albedo(:,:,i_r_diffuse,i_r_ir ), 'SFC_ALB_IR_dif' , 'surface albedo (IR; diffuse)', '1' , fill_halo=.true. )
576  call file_history_in( sfc_albedo(:,:,i_r_direct ,i_r_nir), 'SFC_ALB_NIR_dir', 'surface albedo (NIR; direct', '1' , fill_halo=.true. )
577  call file_history_in( sfc_albedo(:,:,i_r_diffuse,i_r_nir), 'SFC_ALB_NIR_dif', 'surface albedo (NIR; diffuse', '1' , fill_halo=.true. )
578  call file_history_in( sfc_albedo(:,:,i_r_direct ,i_r_vis), 'SFC_ALB_VIS_dir', 'surface albedo (VIS; direct', '1' , fill_halo=.true. )
579  call file_history_in( sfc_albedo(:,:,i_r_diffuse,i_r_vis), 'SFC_ALB_VIS_dif', 'surface albedo (VIS; diffuse', '1' , fill_halo=.true. )
580  call file_history_in( sfc_z0m(:,:), 'SFC_Z0M', 'roughness length (momentum)', 'm' , fill_halo=.true. )
581  call file_history_in( sfc_z0h(:,:), 'SFC_Z0H', 'roughness length (heat)', 'm' , fill_halo=.true. )
582  call file_history_in( sfc_z0e(:,:), 'SFC_Z0E', 'roughness length (vapor)', 'm' , fill_halo=.true. )
583  call file_history_in( sflx_mw(:,:), 'MWFLX', 'w-momentum flux', 'kg/m/s2' )
584  call file_history_in( sflx_mu(:,:), 'MUFLX', 'u-momentum flux', 'kg/m/s2' )
585  call file_history_in( sflx_mv(:,:), 'MVFLX', 'v-momentum flux', 'kg/m/s2' )
586  call file_history_in( sflx_sh(:,:), 'SHFLX', 'sensible heat flux', 'W/m2' , fill_halo=.true. )
587  call file_history_in( sflx_lh(:,:), 'LHFLX', 'latent heat flux', 'W/m2' , fill_halo=.true. )
588  call file_history_in( sflx_gh(:,:), 'GHFLX', 'ground heat flux (downward)', 'W/m2' , fill_halo=.true. )
589  do iq = 1, qa
590  call file_history_in( sflx_qtrc(:,:,iq), 'SFLX_'//trim(tracer_name(iq)), &
591  'surface '//trim(tracer_name(iq))//' flux', &
592  'kg/m2/s' , fill_halo=.true. )
593  enddo
594  call file_history_in( sflx_engi(:,:), 'SFLX_ENGI', 'ground internal energy flux (merged)', 'W/m2' , fill_halo=.true. )
595 
596  call file_history_in( ustar(:,:), 'Ustar', 'friction velocity', 'm/s' , fill_halo=.true. )
597  call file_history_in( tstar(:,:), 'Tstar', 'temperature scale', 'K' , fill_halo=.true. )
598  call file_history_in( qstar(:,:), 'Qstar', 'moisuter scale', 'kg/kg', fill_halo=.true. )
599  call file_history_in( wstar(:,:), 'Wstar', 'convective velocity scale', 'm/s', fill_halo=.true. )
600  call file_history_in( rlmo(:,:), 'RLmo', 'inverse of Obukhov length', '1/m' , fill_halo=.true. )
601 
602  call file_history_in( u10(:,:), 'U10', '10m x-wind', 'm/s' , fill_halo=.true. )
603  call file_history_in( v10(:,:), 'V10', '10m y-wind', 'm/s' , fill_halo=.true. )
604  call file_history_in( t2(:,:), 'T2 ', '2m air temperature', 'K' , fill_halo=.true. )
605  call file_history_in( q2(:,:), 'Q2 ', '2m specific humidity', 'kg/kg', fill_halo=.true. )
606 
607 
608  !$acc data create(work)
609 
610  call file_history_query( hist_uabs10, do_put )
611  if ( do_put ) then
612  !$omp parallel do
613  !$acc kernels
614  do j = js, je
615  do i = is, ie
616  work(i,j) = sqrt( u10(i,j)**2 + v10(i,j)**2 )
617  enddo
618  enddo
619  !$acc end kernels
620  call file_history_put( hist_uabs10, work(:,:) )
621  end if
622 
623  call file_history_query( hist_u10m, do_put )
624  if ( do_put ) then
625  !$omp parallel do
626  !$acc kernels
627  do j = js, je
628  do i = is, ie
629  work(i,j) = u10(i,j) * rotc(i,j,1) - v10(i,j) * rotc(i,j,2)
630  enddo
631  enddo
632  !$acc end kernels
633  call file_history_put( hist_u10m, work(:,:) )
634  end if
635 
636  call file_history_query( hist_v10m, do_put )
637  if ( do_put ) then
638  !$omp parallel do
639  !$acc kernels
640  do j = js, je
641  do i = is, ie
642  work(i,j) = u10(i,j) * rotc(i,j,2) + v10(i,j) * rotc(i,j,1)
643  enddo
644  enddo
645  !$acc end kernels
646  call file_history_put( hist_v10m, work(:,:) )
647  end if
648 
649  call file_history_query( hist_rh2, do_put )
650  if ( do_put ) then
651  call atmos_saturation_psat_liq( &
652  ia, is, ie, ja, js, je, &
653  t2(:,:), & ! (in)
654  work(:,:) ) ! (out)
655  !$omp parallel do
656  !$acc kernels
657  do j = js, je
658  do i = is, ie
659  work(i,j) = sfc_dens(i,j) * q2(i,j) &
660  / work(i,j) * rvap * t2(i,j) &
661  * 100.0_rp
662  enddo
663  enddo
664  !$acc end kernels
665  call file_history_put( hist_rh2, work(:,:) )
666  end if
667 
668 
669  call file_history_query( hist_mslp, do_put )
670  if ( do_put ) then
671  call barometric_law_mslp( ka, ks, ke, ia, is, ie, ja, js, je, & ! [IN]
672  pres(:,:,:), temp(:,:,:), qv(:,:,:), & ! [IN]
673  real_cz(:,:,:), & ! [IN]
674  work(:,:) ) ! [OUT]
675  call file_history_put( hist_mslp, work(:,:) )
676  end if
677 
678  !$acc end data
679 
680  return
681  end subroutine history_output
682 
683 end module mod_atmos_phy_sf_driver
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
Definition: mod_atmos_phy_sf_vars.F90:67
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_dens
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
Definition: mod_atmos_phy_sf_vars.F90:72
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
Definition: mod_atmos_phy_sf_vars.F90:78
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_cpl_sfc_index::i_r_direct
integer, parameter, public i_r_direct
Definition: scale_cpl_sfc_index.F90:37
scale_const::const_epstvap
real(rp), public const_epstvap
1 / epsilon - 1
Definition: scale_const.F90:76
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_flux
subroutine, public atmos_phy_sf_bulk_flux(IA, IS, IE, JA, JS, JE, ATM_W, ATM_U, ATM_V, ATM_TEMP, ATM_PRES, ATM_QV, SFC_DENS, SFC_TEMP, SFC_PRES, SFC_Z0M, SFC_Z0H, SFC_Z0E, PBL, ATM_Z1, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QV, Ustar, Tstar, Qstar, Wstar, RLmo, U10, V10, T2, Q2)
Calculate surface flux.
Definition: scale_atmos_phy_sf_bulk.F90:95
mod_atmos_vars::rhoq_tp
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
Definition: mod_atmos_vars.F90:121
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:39
mod_atmos_phy_sf_vars::atmos_phy_sf_v10
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_v10
Definition: mod_atmos_phy_sf_vars.F90:96
mod_atmos_vars::pott
real(rp), dimension(:,:,:), allocatable, target, public pott
Definition: mod_atmos_vars.F90:133
scale_tracer::tracer_mass
real(rp), dimension(qa_max), public tracer_mass
Definition: scale_tracer.F90:47
mod_atmos_phy_sf_driver::history_output
subroutine history_output
Definition: mod_atmos_phy_sf_driver.F90:516
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0e
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
Definition: mod_atmos_phy_sf_vars.F90:69
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_z1
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_z1
Height of the lowermost grid from surface (cell center) [m].
Definition: scale_atmos_grid_cartesC_real.F90:60
scale_cpl_sfc_index::i_r_diffuse
integer, parameter, public i_r_diffuse
Definition: scale_cpl_sfc_index.F90:38
scale_topography::topography_tansl_y
real(rp), dimension(:,:), allocatable, public topography_tansl_y
tan(slope_y)
Definition: scale_topography.F90:41
scale_precision
module PRECISION
Definition: scale_precision.F90:14
mod_atmos_vars::rhov_tp
real(rp), dimension(:,:,:), allocatable, public rhov_tp
Definition: mod_atmos_vars.F90:118
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
mod_atmos_phy_bl_vars::atmos_phy_bl_zi
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_zi
Definition: mod_atmos_phy_bl_vars.F90:67
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_qstar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_qstar
Definition: mod_atmos_phy_sf_vars.F90:92
mod_atmos_phy_sf_vars::atmos_phy_sf_rhoh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhoh
Definition: mod_atmos_phy_sf_vars.F90:62
mod_atmos_phy_sf_vars::atmos_phy_sf_ustar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_ustar
Definition: mod_atmos_phy_sf_vars.F90:90
scale_const::const_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:68
mod_atmos_phy_sf_vars::atmos_phy_sf_wstar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_wstar
Definition: mod_atmos_phy_sf_vars.F90:93
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totarea
real(rp), public atmos_grid_cartesc_real_totarea
total area (xy, local) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:78
mod_atmos_vars::cptot
real(rp), dimension(:,:,:), allocatable, target, public cptot
Definition: mod_atmos_vars.F90:143
scale_cpl_sfc_index::i_r_ir
integer, parameter, public i_r_ir
Definition: scale_cpl_sfc_index.F90:29
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
scale_bulkflux
module Surface bulk flux
Definition: scale_bulkflux.F90:12
mod_atmos_vars::rhot_av
real(rp), dimension(:,:,:), pointer, public rhot_av
Definition: mod_atmos_vars.F90:94
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
Definition: mod_atmos_phy_sf_vars.F90:79
mod_atmos_phy_rd_vars
module Atmosphere / Physics Radiation
Definition: mod_atmos_phy_rd_vars.F90:12
mod_atmos_phy_sf_vars
module ATMOSPHERIC Surface Variables
Definition: mod_atmos_phy_sf_vars.F90:12
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
Definition: mod_atmos_phy_sf_vars.F90:77
scale_topography
module TOPOGRAPHY
Definition: scale_topography.F90:11
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:114
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:80
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
Definition: mod_atmos_phy_sf_vars.F90:68
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_temp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
Definition: mod_atmos_phy_sf_vars.F90:65
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_tracer::tracer_engi0
real(rp), dimension(qa_max), public tracer_engi0
Definition: scale_tracer.F90:45
scale_atmos_phy_sf_bulk
module atmosphere / physics / surface / bulk
Definition: scale_atmos_phy_sf_bulk.F90:13
scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_setup
subroutine, public atmos_phy_sf_bulk_setup
Setup.
Definition: scale_atmos_phy_sf_bulk.F90:51
scale_const::const_cpvap
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:69
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qvex
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_qvex
Definition: mod_atmos_phy_sf_vars.F90:84
mod_atmos_phy_sf_vars::atmos_phy_sf_u10
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10
Definition: mod_atmos_phy_sf_vars.F90:95
mod_atmos_phy_sf_driver::atmos_phy_sf_driver_calc_tendency
subroutine, public atmos_phy_sf_driver_calc_tendency(update_flag)
calculation tendency
Definition: mod_atmos_phy_sf_driver.F90:149
mod_atmos_vars::rhou_tp
real(rp), dimension(:,:,:), allocatable, public rhou_tp
Definition: mod_atmos_vars.F90:117
mod_atmos_admin::atmos_phy_sf_type
character(len=h_short), public atmos_phy_sf_type
Definition: mod_atmos_admin.F90:40
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_io
module STDIO
Definition: scale_io.F90:10
mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_lw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_dn
Definition: mod_atmos_phy_rd_vars.F90:62
mod_atmos_phy_bl_vars
module atmosphere / physics / PBL
Definition: mod_atmos_phy_bl_vars.F90:12
mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_sw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
Definition: mod_atmos_phy_rd_vars.F90:64
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:76
scale_cpl_sfc_index::i_r_nir
integer, parameter, public i_r_nir
Definition: scale_cpl_sfc_index.F90:30
mod_atmos_phy_sf_vars::atmos_phy_sf_rhoq_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_rhoq_t
Definition: mod_atmos_phy_sf_vars.F90:63
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_shex
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_shex
Definition: mod_atmos_phy_sf_vars.F90:82
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
Definition: mod_atmos_phy_sf_vars.F90:66
mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup
subroutine, public atmos_phy_sf_driver_setup
Setup.
Definition: mod_atmos_phy_sf_driver.F90:54
scale_atmos_hydrometeor::lhv
real(rp), public lhv
latent heat of vaporization for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:144
scale_time::time_dtsec_atmos_phy_sf
real(dp), public time_dtsec_atmos_phy_sf
time interval of physics(surface flux) [sec]
Definition: scale_time.F90:40
scale_topography::topography_tansl_x
real(rp), dimension(:,:), allocatable, public topography_tansl_x
tan(slope_x)
Definition: scale_topography.F90:40
scale_tracer::tracer_cv
real(rp), dimension(qa_max), public tracer_cv
Definition: scale_tracer.F90:42
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
mod_atmos_phy_sf_vars::atmos_phy_sf_rlmo
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rlmo
Definition: mod_atmos_phy_sf_vars.F90:100
mod_atmos_vars::v
real(rp), dimension(:,:,:), allocatable, target, public v
Definition: mod_atmos_vars.F90:131
scale_tracer::tracer_name
character(len=h_short), dimension(qa_max), public tracer_name
Definition: scale_tracer.F90:39
mod_atmos_vars::w
real(rp), dimension(:,:,:), allocatable, target, public w
Definition: mod_atmos_vars.F90:129
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_area
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_area
horizontal area ( xy, normal z) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:66
mod_atmos_vars::momz_tp
real(rp), dimension(:,:,:), allocatable, public momz_tp
Definition: mod_atmos_vars.F90:116
scale_atmos_phy_sf_const::atmos_phy_sf_const_setup
subroutine, public atmos_phy_sf_const_setup
Setup.
Definition: scale_atmos_phy_sf_const.F90:65
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
mod_atmos_vars::temp
real(rp), dimension(:,:,:), allocatable, target, public temp
Definition: mod_atmos_vars.F90:134
mod_atmos_vars::dens_tp
real(rp), dimension(:,:,:), allocatable, public dens_tp
Definition: mod_atmos_vars.F90:115
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_time
module TIME
Definition: scale_time.F90:11
mod_atmos_vars::qv
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
Definition: mod_atmos_vars.F90:97
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, target, public atmos_phy_sf_sflx_qtrc
Definition: mod_atmos_phy_sf_vars.F90:86
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
mod_atmos_vars::pres
real(rp), dimension(:,:,:), allocatable, target, public pres
Definition: mod_atmos_vars.F90:135
mod_atmos_phy_sf_driver
module ATMOSPHERE / Physics Surface fluxes
Definition: mod_atmos_phy_sf_driver.F90:12
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_engi
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_engi
Definition: mod_atmos_phy_sf_vars.F90:87
mod_atmos_vars::dens_av
real(rp), dimension(:,:,:), pointer, public dens_av
Definition: mod_atmos_vars.F90:90
mod_atmos_vars::u
real(rp), dimension(:,:,:), allocatable, target, public u
Definition: mod_atmos_vars.F90:130
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
mod_atmos_phy_sf_vars::atmos_phy_sf_t2
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
Definition: mod_atmos_phy_sf_vars.F90:97
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_lh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
Definition: mod_atmos_phy_sf_vars.F90:81
scale_tracer::tracer_cp
real(rp), dimension(qa_max), public tracer_cp
Definition: scale_tracer.F90:43
mod_atmos_phy_sf_vars::atmos_phy_sf_dens_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_dens_t
Definition: mod_atmos_phy_sf_vars.F90:58
scale_atmos_hydrostatic
module atmosphere / hydrostatic barance
Definition: scale_atmos_hydrostatic.F90:12
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_pres
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
Definition: mod_atmos_phy_sf_vars.F90:73
mod_atmos_vars::rhoh_p
real(rp), dimension(:,:,:), allocatable, public rhoh_p
Definition: mod_atmos_vars.F90:120
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:109
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
scale_cpl_sfc_index
module coupler / surface-atmospehre
Definition: scale_cpl_sfc_index.F90:11
scale_atmos_bottom
module atmosphere / bottom boundary extrapolation
Definition: scale_atmos_bottom.F90:12
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_cpl_sfc_index::i_r_vis
integer, parameter, public i_r_vis
Definition: scale_cpl_sfc_index.F90:31
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
Definition: mod_atmos_phy_sf_vars.F90:80
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
scale_file_history::file_history_reg
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
Definition: scale_file_history.F90:685
mod_atmos_phy_sf_vars::atmos_phy_sf_q2
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
Definition: mod_atmos_phy_sf_vars.F90:98
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:43
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_gh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh
Definition: mod_atmos_phy_sf_vars.F90:85
mod_atmos_phy_sf_vars::atmos_phy_sf_rhou_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhou_t
Definition: mod_atmos_phy_sf_vars.F90:60
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_atmos_phy_sf_const::atmos_phy_sf_const_flux
subroutine, public atmos_phy_sf_const_flux(IA, IS, IE, JA, JS, JE, ATM_W, ATM_U, ATM_V, ATM_TEMP, ATM_Z1, SFC_DENS, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QV, U10, V10)
Constant flux.
Definition: scale_atmos_phy_sf_const.F90:110
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
mod_cpl_admin::cpl_sw
logical, public cpl_sw
Definition: mod_cpl_admin.F90:33
mod_atmos_vars::cvtot
real(rp), dimension(:,:,:), allocatable, target, public cvtot
Definition: mod_atmos_vars.F90:142
scale_atmos_bottom::atmos_bottom_estimate
subroutine, public atmos_bottom_estimate(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, QV, SFC_TEMP, FZ, SFC_DENS, SFC_PRES)
Calc bottom boundary of atmosphere (just above surface)
Definition: scale_atmos_bottom.F90:51
scale_atmos_phy_sf_const
module atmosphere / physics / surface / const
Definition: scale_atmos_phy_sf_const.F90:13
mod_atmos_phy_sf_vars::atmos_phy_sf_tstar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_tstar
Definition: mod_atmos_phy_sf_vars.F90:91
mod_cpl_admin
module Coupler admin
Definition: mod_cpl_admin.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_momz_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momz_t
Definition: mod_atmos_phy_sf_vars.F90:59
mod_atmos_phy_sf_vars::atmos_phy_sf_rhov_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhov_t
Definition: mod_atmos_phy_sf_vars.F90:61
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
Definition: scale_atmos_grid_cartesC_metric.F90:36
mod_atmos_admin::atmos_sw_phy_sf
logical, public atmos_sw_phy_sf
Definition: mod_atmos_admin.F90:56