SCALE-RM
mod_atmos_phy_sf_vars.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
16  !-----------------------------------------------------------------------------
17  !
18  !++ used modules
19  !
20  use scale_precision
21  use scale_stdio
22  use scale_prof
24  use scale_tracer
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public procedure
31  !
32  public :: atmos_phy_sf_vars_setup
37 
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public parameters & variables
41  !
42  logical, public :: atmos_phy_sf_restart_output = .false.
43 
44  character(len=H_LONG), public :: atmos_phy_sf_restart_in_basename = ''
45  character(len=H_LONG), public :: atmos_phy_sf_restart_out_basename = ''
46  character(len=H_MID), public :: atmos_phy_sf_restart_out_title = 'ATMOS_PHY_SF restart'
47  character(len=H_MID), public :: atmos_phy_sf_restart_out_dtype = 'DEFAULT'
48 
49  real(RP), public, allocatable :: atmos_phy_sf_dens_t(:,:) ! tendency DENS [ kg/m3/s]
50  real(RP), public, allocatable :: atmos_phy_sf_momz_t(:,:) ! tendency MOMZ [m/s*kg/m3/s]
51  real(RP), public, allocatable :: atmos_phy_sf_momx_t(:,:) ! tendency MOMX [m/s*kg/m3/s]
52  real(RP), public, allocatable :: atmos_phy_sf_momy_t(:,:) ! tendency MOMY [m/s*kg/m3/s]
53  real(RP), public, allocatable :: atmos_phy_sf_rhot_t(:,:) ! tendency RHOT [K *kg/m3/s]
54  real(RP), public, allocatable :: atmos_phy_sf_rhoq_t(:,:,:) ! tendency rho*QTRC [ kg/kg/s]
55 
56  real(RP), public, allocatable :: atmos_phy_sf_sfc_temp (:,:) ! surface skin temperature [K]
57  real(RP), public, allocatable :: atmos_phy_sf_sfc_albedo(:,:,:) ! surface albedo [0-1]
58  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0m (:,:) ! surface roughness length, ocean only [m]
59  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0h (:,:) ! surface roughness length, ocean only [m]
60  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0e (:,:) ! surface roughness length, ocean only [m]
61 
62  ! surface diagnostic variables
63  real(RP), public, allocatable :: atmos_phy_sf_sfc_dens (:,:) ! surface atmosphere density [kg/m3]
64  real(RP), public, allocatable :: atmos_phy_sf_sfc_pres (:,:) ! surface atmosphere pressure [Pa]
65 
66  real(RP), public, allocatable :: atmos_phy_sf_sflx_mw (:,:) ! z-momentum flux (area center) [m/s*kg/m2/s]
67  real(RP), public, allocatable :: atmos_phy_sf_sflx_mu (:,:) ! x-momentum flux (area center) [m/s*kg/m2/s]
68  real(RP), public, allocatable :: atmos_phy_sf_sflx_mv (:,:) ! y-momentum flux (area center) [m/s*kg/m2/s]
69  real(RP), public, allocatable :: atmos_phy_sf_sflx_sh (:,:) ! sensible heat flux [J/m2/s]
70  real(RP), public, allocatable :: atmos_phy_sf_sflx_lh (:,:) ! latent heat flux [J/m2/s]
71  real(RP), public, allocatable :: atmos_phy_sf_sflx_gh (:,:) ! ground heat flux [J/m2/s]
72  real(RP), public, allocatable :: atmos_phy_sf_sflx_qtrc (:,:,:) ! tracer mass flux [kg/m2/s]
73 
74  real(RP), public, allocatable :: atmos_phy_sf_u10 (:,:) ! 10m x-wind [m/s]
75  real(RP), public, allocatable :: atmos_phy_sf_v10 (:,:) ! 10m y-wind [m/s]
76  real(RP), public, allocatable :: atmos_phy_sf_t2 (:,:) ! 2m temperature [K]
77  real(RP), public, allocatable :: atmos_phy_sf_q2 (:,:) ! 2m specific humidity [kg/kg]
78 
79 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_QEMIS(:,:,:) ! tracer emission flux [kg/m2/s]
80 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_QDEP (:,:,:) ! tracer deposition flux [kg/m2/s]
81 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_VDEP (:,:,:) ! tracer deposition velocity [m/s]
82 
83  !-----------------------------------------------------------------------------
84  !
85  !++ Private procedure
86  !
87  !-----------------------------------------------------------------------------
88  !
89  !++ Private parameters & variables
90  !
91  integer, private, parameter :: vmax = 6
92  integer, private, parameter :: i_sfc_temp = 1
93  integer, private, parameter :: i_sfc_albedo_lw = 2
94  integer, private, parameter :: i_sfc_albedo_sw = 3
95  integer, private, parameter :: i_sfc_z0m = 4
96  integer, private, parameter :: i_sfc_z0h = 5
97  integer, private, parameter :: i_sfc_z0e = 6
98 
99  character(len=H_SHORT), private :: var_name(vmax)
100  character(len=H_MID), private :: var_desc(vmax)
101  character(len=H_SHORT), private :: var_unit(vmax)
102 
103  data var_name / 'SFC_TEMP', &
104  'SFC_ALB_LW', &
105  'SFC_ALB_SW', &
106  'SFC_Z0M', &
107  'SFC_Z0H', &
108  'SFC_Z0E' /
109 
110  data var_desc / 'surface skin temperature', &
111  'surface albedo for longwave', &
112  'surface albedo for shortwave', &
113  'surface roughness length (momentum)', &
114  'surface roughness length (heat)', &
115  'surface roughness length (vapor)' /
116 
117  data var_unit / 'K', &
118  '0-1', &
119  '0-1', &
120  'm', &
121  'm', &
122  'm' /
123 
124  real(RP), private :: atmos_phy_sf_default_sfc_temp = 300.0_rp
125  real(RP), private :: atmos_phy_sf_default_sfc_albedo = 0.4_rp
126  real(RP), private :: atmos_phy_sf_default_sfc_z0 = 1e-4_rp
127 
128  !-----------------------------------------------------------------------------
129 contains
130  !-----------------------------------------------------------------------------
132  subroutine atmos_phy_sf_vars_setup
133  use scale_process, only: &
135  use scale_const, only: &
136  undef => const_undef
137  implicit none
138 
139  namelist / param_atmos_phy_sf_vars / &
145  atmos_phy_sf_default_sfc_temp, &
146  atmos_phy_sf_default_sfc_albedo, &
147  atmos_phy_sf_default_sfc_z0
148 
149  integer :: ierr
150  integer :: iv
151  !---------------------------------------------------------------------------
152 
153  if( io_l ) write(io_fid_log,*)
154  if( io_l ) write(io_fid_log,*) '++++++ Module[VARS] / Categ[ATMOS PHY_SF] / Origin[SCALE-RM]'
155 
156  allocate( atmos_phy_sf_dens_t(ia,ja) )
157  allocate( atmos_phy_sf_momz_t(ia,ja) )
158  allocate( atmos_phy_sf_momx_t(ia,ja) )
159  allocate( atmos_phy_sf_momy_t(ia,ja) )
160  allocate( atmos_phy_sf_rhot_t(ia,ja) )
161  allocate( atmos_phy_sf_rhoq_t(ia,ja,qa) )
162  atmos_phy_sf_dens_t(:,:) = undef
163  atmos_phy_sf_momz_t(:,:) = undef
164  atmos_phy_sf_momx_t(:,:) = undef
165  atmos_phy_sf_momy_t(:,:) = undef
166  atmos_phy_sf_rhot_t(:,:) = undef
167  atmos_phy_sf_rhoq_t(:,:,:) = undef
168 
169  allocate( atmos_phy_sf_sfc_temp(ia,ja) )
170  allocate( atmos_phy_sf_sfc_albedo(ia,ja,2) )
171  allocate( atmos_phy_sf_sfc_z0m(ia,ja) )
172  allocate( atmos_phy_sf_sfc_z0h(ia,ja) )
173  allocate( atmos_phy_sf_sfc_z0e(ia,ja) )
174 ! ATMOS_PHY_SF_SFC_TEMP (:,:) = UNDEF ! [del] 2014/8/28 A.Noda
175 ! ATMOS_PHY_SF_SFC_albedo(:,:,:) = UNDEF ! [del] 2014/8/28 A.Noda
176  atmos_phy_sf_sfc_z0m(:,:) = undef
177  atmos_phy_sf_sfc_z0h(:,:) = undef
178  atmos_phy_sf_sfc_z0e(:,:) = undef
179 
180  allocate( atmos_phy_sf_sfc_dens(ia,ja) )
181  allocate( atmos_phy_sf_sfc_pres(ia,ja) )
182  atmos_phy_sf_sfc_dens(:,:) = undef
183  atmos_phy_sf_sfc_pres(:,:) = undef
184 
185  allocate( atmos_phy_sf_sflx_mw(ia,ja) )
186  allocate( atmos_phy_sf_sflx_mu(ia,ja) )
187  allocate( atmos_phy_sf_sflx_mv(ia,ja) )
188  allocate( atmos_phy_sf_sflx_sh(ia,ja) )
189  allocate( atmos_phy_sf_sflx_lh(ia,ja) )
190  allocate( atmos_phy_sf_sflx_gh(ia,ja) )
191  allocate( atmos_phy_sf_sflx_qtrc(ia,ja,qa) )
192  atmos_phy_sf_sflx_mw(:,:) = undef
193  atmos_phy_sf_sflx_mu(:,:) = undef
194  atmos_phy_sf_sflx_mv(:,:) = undef
195  atmos_phy_sf_sflx_sh(:,:) = undef
196  atmos_phy_sf_sflx_lh(:,:) = undef
197  atmos_phy_sf_sflx_gh(:,:) = undef
198  atmos_phy_sf_sflx_qtrc(:,:,:) = undef
199 
200  allocate( atmos_phy_sf_u10(ia,ja) )
201  allocate( atmos_phy_sf_v10(ia,ja) )
202  allocate( atmos_phy_sf_t2(ia,ja) )
203  allocate( atmos_phy_sf_q2(ia,ja) )
204  atmos_phy_sf_u10(:,:) = undef
205  atmos_phy_sf_v10(:,:) = undef
206  atmos_phy_sf_t2(:,:) = undef
207  atmos_phy_sf_q2(:,:) = undef
208 
209  !--- read namelist
210  rewind(io_fid_conf)
211  read(io_fid_conf,nml=param_atmos_phy_sf_vars,iostat=ierr)
212  if( ierr < 0 ) then !--- missing
213  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
214  elseif( ierr > 0 ) then !--- fatal error
215  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_VARS. Check!'
216  call prc_mpistop
217  endif
218  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_sf_vars)
219 
220  ! [add] 2014/08/28 A.Noda
221  atmos_phy_sf_sfc_temp(:,:) = atmos_phy_sf_default_sfc_temp
222  atmos_phy_sf_sfc_albedo(:,:,:) = atmos_phy_sf_default_sfc_albedo
223  atmos_phy_sf_sfc_z0m(:,:) = atmos_phy_sf_default_sfc_z0
224  atmos_phy_sf_sfc_z0h(:,:) = atmos_phy_sf_default_sfc_z0
225  atmos_phy_sf_sfc_z0e(:,:) = atmos_phy_sf_default_sfc_z0
226 
227  if( io_l ) write(io_fid_log,*)
228  if( io_l ) write(io_fid_log,*) '*** [ATMOS_PHY_SF] prognostic/diagnostic variables'
229  if( io_l ) write(io_fid_log,'(1x,A,A15,A,A32,3(A))') &
230  '*** |','VARNAME ','|', 'DESCRIPTION ','[', 'UNIT ',']'
231  do iv = 1, vmax
232  if( io_l ) write(io_fid_log,'(1x,A,i3,A,A15,A,A32,3(A))') &
233  '*** NO.',iv,'|',var_name(iv),'|',var_desc(iv),'[',var_unit(iv),']'
234  enddo
235 
236  if( io_l ) write(io_fid_log,*)
237  if ( atmos_phy_sf_restart_in_basename /= '' ) then
238  if( io_l ) write(io_fid_log,*) '*** Restart input? : ', trim(atmos_phy_sf_restart_in_basename)
239  else
240  if( io_l ) write(io_fid_log,*) '*** Restart input? : NO'
241  endif
243  .AND. atmos_phy_sf_restart_out_basename /= '' ) then
244  if( io_l ) write(io_fid_log,*) '*** Restart output? : ', trim(atmos_phy_sf_restart_out_basename)
245  else
246  if( io_l ) write(io_fid_log,*) '*** Restart output? : NO'
248  endif
249 
250  return
251  end subroutine atmos_phy_sf_vars_setup
252 
253  !-----------------------------------------------------------------------------
255  subroutine atmos_phy_sf_vars_fillhalo
256  use scale_const, only: &
257  i_sw => const_i_sw, &
258  i_lw => const_i_lw
259  use scale_comm, only: &
260  comm_vars8, &
261  comm_wait
262  implicit none
263  !---------------------------------------------------------------------------
264 
265  call comm_vars8( atmos_phy_sf_sfc_temp(:,:), 1 )
266  call comm_vars8( atmos_phy_sf_sfc_albedo(:,:,i_lw), 2 )
267  call comm_vars8( atmos_phy_sf_sfc_albedo(:,:,i_sw), 3 )
268  call comm_vars8( atmos_phy_sf_sfc_z0m(:,:), 4 )
269  call comm_vars8( atmos_phy_sf_sfc_z0h(:,:), 5 )
270  call comm_vars8( atmos_phy_sf_sfc_z0e(:,:), 6 )
271 
272  call comm_wait ( atmos_phy_sf_sfc_temp(:,:), 1 )
273  call comm_wait ( atmos_phy_sf_sfc_albedo(:,:,i_lw), 2 )
274  call comm_wait ( atmos_phy_sf_sfc_albedo(:,:,i_sw), 3 )
275  call comm_wait ( atmos_phy_sf_sfc_z0m(:,:), 4 )
276  call comm_wait ( atmos_phy_sf_sfc_z0h(:,:), 5 )
277  call comm_wait ( atmos_phy_sf_sfc_z0e(:,:), 6 )
278 
279  return
280  end subroutine atmos_phy_sf_vars_fillhalo
281 
282  !-----------------------------------------------------------------------------
285  use scale_const, only: &
286  i_sw => const_i_sw, &
287  i_lw => const_i_lw
288  use scale_fileio, only: &
289  fileio_read
290  use scale_rm_statistics, only: &
291  stat_total
292  implicit none
293 
294  real(RP) :: total
295  !---------------------------------------------------------------------------
296 
297  if( io_l ) write(io_fid_log,*)
298  if( io_l ) write(io_fid_log,*) '*** Input restart file (ATMOS_PHY_SF) ***'
299 
300  if ( atmos_phy_sf_restart_in_basename /= '' ) then
301  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(atmos_phy_sf_restart_in_basename)
302 
303  call fileio_read( atmos_phy_sf_sfc_temp(:,:), & ! [OUT]
304  atmos_phy_sf_restart_in_basename, var_name(1), 'XY', step=1 ) ! [IN]
305  call fileio_read( atmos_phy_sf_sfc_albedo(:,:,i_lw), & ! [OUT]
306  atmos_phy_sf_restart_in_basename, var_name(2), 'XY', step=1 ) ! [IN]
307  call fileio_read( atmos_phy_sf_sfc_albedo(:,:,i_sw), & ! [OUT]
308  atmos_phy_sf_restart_in_basename, var_name(3), 'XY', step=1 ) ! [IN]
309  call fileio_read( atmos_phy_sf_sfc_z0m(:,:), & ! [OUT]
310  atmos_phy_sf_restart_in_basename, var_name(4), 'XY', step=1 ) ! [IN]
311  call fileio_read( atmos_phy_sf_sfc_z0h(:,:), & ! [OUT]
312  atmos_phy_sf_restart_in_basename, var_name(5), 'XY', step=1 ) ! [IN]
313  call fileio_read( atmos_phy_sf_sfc_z0e(:,:), & ! [OUT]
314  atmos_phy_sf_restart_in_basename, var_name(6), 'XY', step=1 ) ! [IN]
315 
317 
318  call stat_total( total, atmos_phy_sf_sfc_temp(:,:), var_name(1) )
319  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_lw), var_name(2) )
320  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_sw), var_name(3) )
321  call stat_total( total, atmos_phy_sf_sfc_z0m(:,:), var_name(4) )
322  call stat_total( total, atmos_phy_sf_sfc_z0h(:,:), var_name(5) )
323  call stat_total( total, atmos_phy_sf_sfc_z0e(:,:), var_name(6) )
324  else
325  if( io_l ) write(io_fid_log,*) '*** restart file for ATMOS_PHY_SF is not specified.'
326  atmos_phy_sf_sfc_temp(:,:) = atmos_phy_sf_default_sfc_temp
327  atmos_phy_sf_sfc_albedo(:,:,:) = atmos_phy_sf_default_sfc_albedo
328  atmos_phy_sf_sfc_z0m(:,:) = atmos_phy_sf_default_sfc_z0
329  atmos_phy_sf_sfc_z0h(:,:) = atmos_phy_sf_default_sfc_z0
330  atmos_phy_sf_sfc_z0e(:,:) = atmos_phy_sf_default_sfc_z0
331  endif
332 
333  return
334  end subroutine atmos_phy_sf_vars_restart_read
335 
336  !-----------------------------------------------------------------------------
339  use scale_const, only: &
340  i_sw => const_i_sw, &
341  i_lw => const_i_lw
342  use scale_time, only: &
344  use scale_fileio, only: &
345  fileio_write
346  use scale_rm_statistics, only: &
347  stat_total
348  implicit none
349 
350  character(len=20) :: timelabel
351  character(len=H_LONG) :: basename
352 
353  real(RP) :: total
354  !---------------------------------------------------------------------------
355 
356  if ( atmos_phy_sf_restart_out_basename /= '' ) then
357 
358  call time_gettimelabel( timelabel )
359  write(basename,'(A,A,A)') trim(atmos_phy_sf_restart_out_basename), '_', trim(timelabel)
360 
361  if( io_l ) write(io_fid_log,*)
362  if( io_l ) write(io_fid_log,*) '*** Output restart file (ATMOS_PHY_SF) ***'
363  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(basename)
364 
366 
367  call stat_total( total, atmos_phy_sf_sfc_temp(:,:), var_name(1) )
368  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_lw), var_name(2) )
369  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_sw), var_name(3) )
370  call stat_total( total, atmos_phy_sf_sfc_z0m(:,:), var_name(4) )
371  call stat_total( total, atmos_phy_sf_sfc_z0h(:,:), var_name(5) )
372  call stat_total( total, atmos_phy_sf_sfc_z0e(:,:), var_name(6) )
373 
374  call fileio_write( atmos_phy_sf_sfc_temp(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
375  var_name(1), var_desc(1), var_unit(1), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
376  call fileio_write( atmos_phy_sf_sfc_albedo(:,:,i_lw), basename, atmos_phy_sf_restart_out_title, & ! [IN]
377  var_name(2), var_desc(2), var_unit(2), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
378  call fileio_write( atmos_phy_sf_sfc_albedo(:,:,i_sw), basename, atmos_phy_sf_restart_out_title, & ! [IN]
379  var_name(3), var_desc(3), var_unit(3), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
380  call fileio_write( atmos_phy_sf_sfc_z0m(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
381  var_name(4), var_desc(4), var_unit(4), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
382  call fileio_write( atmos_phy_sf_sfc_z0h(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
383  var_name(5), var_desc(5), var_unit(5), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
384  call fileio_write( atmos_phy_sf_sfc_z0e(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
385  var_name(6), var_desc(6), var_unit(6), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
386  endif
387 
388  return
389  end subroutine atmos_phy_sf_vars_restart_write
390 
391  !-----------------------------------------------------------------------------
393  subroutine atmos_phy_sf_vars_external_in( &
394  sfc_temp_in, &
395  albedo_lw_in, &
396  albedo_sw_in, &
397  sfc_z0m_in, &
398  sfc_z0h_in, &
399  sfc_z0e_in )
400  use scale_const, only: &
401  i_sw => const_i_sw, &
402  i_lw => const_i_lw
403  implicit none
404 
405  real(RP), intent(in) :: sfc_temp_in (:,:)
406  real(RP), intent(in) :: albedo_lw_in(:,:)
407  real(RP), intent(in) :: albedo_sw_in(:,:)
408  real(RP), intent(in) :: sfc_z0m_in (:,:)
409  real(RP), intent(in) :: sfc_z0h_in (:,:)
410  real(RP), intent(in) :: sfc_z0e_in (:,:)
411  !---------------------------------------------------------------------------
412 
413  if( io_l ) write(io_fid_log,*)
414  if( io_l ) write(io_fid_log,*) '*** External Input (PHY_SF) ***'
415 
416  atmos_phy_sf_sfc_temp(:,:) = sfc_temp_in(:,:)
417  atmos_phy_sf_sfc_albedo(:,:,i_lw) = albedo_lw_in(:,:)
418  atmos_phy_sf_sfc_albedo(:,:,i_sw) = albedo_sw_in(:,:)
419  atmos_phy_sf_sfc_z0m(:,:) = sfc_z0m_in(:,:)
420  atmos_phy_sf_sfc_z0h(:,:) = sfc_z0h_in(:,:)
421  atmos_phy_sf_sfc_z0e(:,:) = sfc_z0e_in(:,:)
422 
424 
425  return
426  end subroutine atmos_phy_sf_vars_external_in
427 
428 end module mod_atmos_phy_sf_vars
subroutine, public atmos_phy_sf_vars_restart_read
Read restart.
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 atmos_phy_sf_vars_external_in(sfc_temp_in, albedo_lw_in, albedo_sw_in, sfc_z0m_in, sfc_z0h_in, sfc_z0e_in)
Input from External I/O.
subroutine, public prc_mpistop
Abort MPI.
character(len=h_long), public atmos_phy_sf_restart_out_basename
basename of the output file
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_rhoq_t
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 STDIO
Definition: scale_stdio.F90:12
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhot_t
integer, public qa
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
module FILE I/O (netcdf)
subroutine, public atmos_phy_sf_vars_setup
Setup.
character(len=h_mid), public atmos_phy_sf_restart_out_dtype
REAL4 or REAL8.
real(rp), public const_undef
Definition: scale_const.F90:43
module Statistics
module grid index
module ATMOSPHERIC Surface Variables
subroutine, public atmos_phy_sf_vars_restart_write
Write restart.
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 atmos_phy_sf_momx_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_dens_t
subroutine, public atmos_phy_sf_vars_fillhalo
HALO Communication.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
subroutine, public time_gettimelabel(timelabel)
generate time label
Definition: scale_time.F90:90
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
character(len=h_long), public atmos_phy_sf_restart_in_basename
basename of the restart file
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
module COMMUNICATION
Definition: scale_comm.F90:23
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
module TIME
Definition: scale_time.F90:15
character(len=h_mid), public atmos_phy_sf_restart_out_title
title of the output file
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momy_t
module profiler
Definition: scale_prof.F90:10
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
module PRECISION
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
integer, public ja
of y whole cells (local, with HALO)
logical, public atmos_phy_sf_restart_output
output restart file?