SCALE-RM
Functions/Subroutines | Variables
mod_atmos_phy_sf_vars Module Reference

module ATMOSPHERIC Surface Variables More...

Functions/Subroutines

subroutine, public atmos_phy_sf_vars_setup
 Setup. More...
 
subroutine, public atmos_phy_sf_vars_fillhalo
 HALO Communication. More...
 
subroutine, public atmos_phy_sf_vars_restart_read
 Read restart. More...
 
subroutine, public atmos_phy_sf_vars_restart_write
 Write restart. More...
 
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. More...
 

Variables

logical, public atmos_phy_sf_restart_output = .false.
 output restart file? More...
 
character(len=h_long), public atmos_phy_sf_restart_in_basename = ''
 basename of the restart file More...
 
character(len=h_long), public atmos_phy_sf_restart_out_basename = ''
 basename of the output file More...
 
character(len=h_mid), public atmos_phy_sf_restart_out_title = 'ATMOS_PHY_SF restart'
 title of the output file More...
 
character(len=h_mid), public atmos_phy_sf_restart_out_dtype = 'DEFAULT'
 REAL4 or REAL8. More...
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_dens_t
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momz_t
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momx_t
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momy_t
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhot_t
 
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_rhoq_t
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
 
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh
 
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sflx_qtrc
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_v10
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
 
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
 

Detailed Description

module ATMOSPHERIC Surface Variables

Description
Container for mod_atmos_phy_sf
Author
Team SCALE
History
  • 2012-03-27 (H.Yashiro) [new]
  • 2013-08-31 (T.Yamaura) [mod]
  • 2014-08-28 (A.Noda) [mod] input of ATMOS_PHY_SF_SFC_{TEMP,albedo}
NAMELIST
  • PARAM_ATMOS_PHY_SF_VARS
    nametypedefault valuecomment
    ATMOS_PHY_SF_RESTART_IN_BASENAME character(len=H_LONG) '' basename of the restart file
    ATMOS_PHY_SF_RESTART_OUTPUT logical .false. output restart file?
    ATMOS_PHY_SF_RESTART_OUT_BASENAME character(len=H_LONG) '' basename of the output file
    ATMOS_PHY_SF_RESTART_OUT_TITLE character(len=H_MID) 'ATMOS_PHY_SF restart' title of the output file
    ATMOS_PHY_SF_RESTART_OUT_DTYPE character(len=H_MID) 'DEFAULT' REAL4 or REAL8
    ATMOS_PHY_SF_DEFAULT_SFC_TEMP real(RP) 300.0_RP
    ATMOS_PHY_SF_DEFAULT_SFC_ALBEDO real(RP) 0.4_RP
    ATMOS_PHY_SF_DEFAULT_SFC_Z0 real(RP) 1E-4_RP

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_sf_vars_setup()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_setup ( )

Setup.

Definition at line 133 of file mod_atmos_phy_sf_vars.f90.

References atmos_phy_sf_dens_t, atmos_phy_sf_momx_t, atmos_phy_sf_momy_t, atmos_phy_sf_momz_t, atmos_phy_sf_q2, atmos_phy_sf_restart_in_basename, atmos_phy_sf_restart_out_basename, atmos_phy_sf_restart_out_dtype, atmos_phy_sf_restart_out_title, atmos_phy_sf_restart_output, atmos_phy_sf_rhoq_t, atmos_phy_sf_rhot_t, atmos_phy_sf_sfc_albedo, atmos_phy_sf_sfc_dens, atmos_phy_sf_sfc_pres, atmos_phy_sf_sfc_temp, atmos_phy_sf_sfc_z0e, atmos_phy_sf_sfc_z0h, atmos_phy_sf_sfc_z0m, atmos_phy_sf_sflx_gh, atmos_phy_sf_sflx_lh, atmos_phy_sf_sflx_mu, atmos_phy_sf_sflx_mv, atmos_phy_sf_sflx_mw, atmos_phy_sf_sflx_qtrc, atmos_phy_sf_sflx_sh, atmos_phy_sf_t2, atmos_phy_sf_u10, atmos_phy_sf_v10, scale_const::const_undef, scale_grid_index::ia, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::ja, scale_process::prc_mpistop(), and scale_tracer::qa.

Referenced by mod_atmos_vars::atmos_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 / &
140  atmos_phy_sf_restart_in_basename, &
141  atmos_phy_sf_restart_output, &
142  atmos_phy_sf_restart_out_basename, &
143  atmos_phy_sf_restart_out_title, &
144  atmos_phy_sf_restart_out_dtype, &
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
242  if ( atmos_phy_sf_restart_output &
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'
247  atmos_phy_sf_restart_output = .false.
248  endif
249 
250  return
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
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, public atmos_phy_sf_v10
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
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
real(rp), public const_undef
Definition: scale_const.F90:43
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
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
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
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
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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_fillhalo()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_fillhalo ( )

HALO Communication.

Definition at line 256 of file mod_atmos_phy_sf_vars.f90.

References atmos_phy_sf_sfc_albedo, atmos_phy_sf_sfc_temp, atmos_phy_sf_sfc_z0e, atmos_phy_sf_sfc_z0h, atmos_phy_sf_sfc_z0m, scale_const::const_i_lw, and scale_const::const_i_sw.

Referenced by atmos_phy_sf_vars_external_in(), atmos_phy_sf_vars_restart_read(), and atmos_phy_sf_vars_restart_write().

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
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
integer, parameter, public i_lw
integer, parameter, public i_sw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
module COMMUNICATION
Definition: scale_comm.F90:23
module CONSTANT
Definition: scale_const.F90:14
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_restart_read()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_read ( )

Read restart.

Definition at line 285 of file mod_atmos_phy_sf_vars.f90.

References atmos_phy_sf_restart_in_basename, atmos_phy_sf_sfc_albedo, atmos_phy_sf_sfc_temp, atmos_phy_sf_sfc_z0e, atmos_phy_sf_sfc_z0h, atmos_phy_sf_sfc_z0m, atmos_phy_sf_vars_fillhalo(), scale_const::const_i_lw, scale_const::const_i_sw, scale_stdio::io_fid_log, and scale_stdio::io_l.

Referenced by mod_atmos_vars::atmos_vars_restart_read().

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 
316  call atmos_phy_sf_vars_fillhalo
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
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
module FILE I/O (netcdf)
integer, parameter, public i_lw
module Statistics
integer, parameter, public i_sw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
module CONSTANT
Definition: scale_const.F90:14
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_restart_write()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_write ( )

Write restart.

Definition at line 339 of file mod_atmos_phy_sf_vars.f90.

References atmos_phy_sf_restart_out_basename, atmos_phy_sf_restart_out_dtype, atmos_phy_sf_restart_out_title, atmos_phy_sf_sfc_albedo, atmos_phy_sf_sfc_temp, atmos_phy_sf_sfc_z0e, atmos_phy_sf_sfc_z0h, atmos_phy_sf_sfc_z0m, atmos_phy_sf_vars_fillhalo(), scale_const::const_i_lw, scale_const::const_i_sw, scale_stdio::io_fid_log, scale_stdio::io_l, and scale_time::time_gettimelabel().

Referenced by mod_atmos_vars::atmos_vars_restart_write().

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 
365  call atmos_phy_sf_vars_fillhalo
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
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
module FILE I/O (netcdf)
integer, parameter, public i_lw
module Statistics
integer, parameter, public i_sw
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
module TIME
Definition: scale_time.F90:15
module CONSTANT
Definition: scale_const.F90:14
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_external_in()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_external_in ( real(rp), dimension (:,:), intent(in)  sfc_temp_in,
real(rp), dimension(:,:), intent(in)  albedo_lw_in,
real(rp), dimension(:,:), intent(in)  albedo_sw_in,
real(rp), dimension (:,:), intent(in)  sfc_z0m_in,
real(rp), dimension (:,:), intent(in)  sfc_z0h_in,
real(rp), dimension (:,:), intent(in)  sfc_z0e_in 
)

Input from External I/O.

Definition at line 400 of file mod_atmos_phy_sf_vars.f90.

References atmos_phy_sf_sfc_albedo, atmos_phy_sf_sfc_temp, atmos_phy_sf_sfc_z0e, atmos_phy_sf_sfc_z0h, atmos_phy_sf_sfc_z0m, atmos_phy_sf_vars_fillhalo(), scale_const::const_i_lw, scale_const::const_i_sw, scale_stdio::io_fid_log, and scale_stdio::io_l.

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 
423  call atmos_phy_sf_vars_fillhalo
424 
425  return
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
integer, parameter, public i_lw
integer, parameter, public i_sw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
module CONSTANT
Definition: scale_const.F90:14
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:

Variable Documentation

◆ atmos_phy_sf_restart_output

logical, public mod_atmos_phy_sf_vars::atmos_phy_sf_restart_output = .false.

output restart file?

Definition at line 42 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_admin_restart::admin_restart_setup(), and atmos_phy_sf_vars_setup().

42  logical, public :: atmos_phy_sf_restart_output = .false.

◆ atmos_phy_sf_restart_in_basename

character(len=h_long), public mod_atmos_phy_sf_vars::atmos_phy_sf_restart_in_basename = ''

basename of the restart file

Definition at line 44 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_admin_restart::admin_restart_setup(), atmos_phy_sf_vars_restart_read(), and atmos_phy_sf_vars_setup().

44  character(len=H_LONG), public :: atmos_phy_sf_restart_in_basename = ''

◆ atmos_phy_sf_restart_out_basename

character(len=h_long), public mod_atmos_phy_sf_vars::atmos_phy_sf_restart_out_basename = ''

basename of the output file

Definition at line 45 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_admin_restart::admin_restart_setup(), atmos_phy_sf_vars_restart_write(), and atmos_phy_sf_vars_setup().

45  character(len=H_LONG), public :: atmos_phy_sf_restart_out_basename = ''

◆ atmos_phy_sf_restart_out_title

character(len=h_mid), public mod_atmos_phy_sf_vars::atmos_phy_sf_restart_out_title = 'ATMOS_PHY_SF restart'

title of the output file

Definition at line 46 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_admin_restart::admin_restart_setup(), atmos_phy_sf_vars_restart_write(), and atmos_phy_sf_vars_setup().

46  character(len=H_MID), public :: atmos_phy_sf_restart_out_title = 'ATMOS_PHY_SF restart'

◆ atmos_phy_sf_restart_out_dtype

character(len=h_mid), public mod_atmos_phy_sf_vars::atmos_phy_sf_restart_out_dtype = 'DEFAULT'

REAL4 or REAL8.

Definition at line 47 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_admin_restart::admin_restart_setup(), atmos_phy_sf_vars_restart_write(), and atmos_phy_sf_vars_setup().

47  character(len=H_MID), public :: atmos_phy_sf_restart_out_dtype = 'DEFAULT'

◆ atmos_phy_sf_dens_t

real(rp), dimension(:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_dens_t

Definition at line 49 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

49  real(RP), public, allocatable :: atmos_phy_sf_dens_t(:,:) ! tendency DENS [ kg/m3/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_dens_t

◆ atmos_phy_sf_momz_t

real(rp), dimension(:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_momz_t

Definition at line 50 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

50  real(RP), public, allocatable :: atmos_phy_sf_momz_t(:,:) ! tendency MOMZ [m/s*kg/m3/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momz_t

◆ atmos_phy_sf_momx_t

real(rp), dimension(:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_momx_t

Definition at line 51 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

51  real(RP), public, allocatable :: atmos_phy_sf_momx_t(:,:) ! tendency MOMX [m/s*kg/m3/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momx_t

◆ atmos_phy_sf_momy_t

real(rp), dimension(:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_momy_t

Definition at line 52 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

52  real(RP), public, allocatable :: atmos_phy_sf_momy_t(:,:) ! tendency MOMY [m/s*kg/m3/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momy_t

◆ atmos_phy_sf_rhot_t

real(rp), dimension(:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_rhot_t

Definition at line 53 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

53  real(RP), public, allocatable :: atmos_phy_sf_rhot_t(:,:) ! tendency RHOT [K *kg/m3/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhot_t

◆ atmos_phy_sf_rhoq_t

real(rp), dimension(:,:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_rhoq_t

Definition at line 54 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

54  real(RP), public, allocatable :: atmos_phy_sf_rhoq_t(:,:,:) ! tendency rho*QTRC [ kg/kg/s]
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_rhoq_t

◆ atmos_phy_sf_sfc_temp

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_temp

◆ atmos_phy_sf_sfc_albedo

real(rp), dimension(:,:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_albedo

◆ atmos_phy_sf_sfc_z0m

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0m

◆ atmos_phy_sf_sfc_z0h

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0h

◆ atmos_phy_sf_sfc_z0e

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0e

◆ atmos_phy_sf_sfc_dens

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_dens

Definition at line 63 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

63  real(RP), public, allocatable :: atmos_phy_sf_sfc_dens (:,:) ! surface atmosphere density [kg/m3]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens

◆ atmos_phy_sf_sfc_pres

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_pres

Definition at line 64 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), and atmos_phy_sf_vars_setup().

64  real(RP), public, allocatable :: atmos_phy_sf_sfc_pres (:,:) ! surface atmosphere pressure [Pa]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres

◆ atmos_phy_sf_sflx_mw

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mw

Definition at line 66 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup(), atmos_phy_sf_vars_setup(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver(), and mod_atmos_driver::atmos_surface_get().

66  real(RP), public, allocatable :: atmos_phy_sf_sflx_mw (:,:) ! z-momentum flux (area center) [m/s*kg/m2/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw

◆ atmos_phy_sf_sflx_mu

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu

Definition at line 67 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup(), atmos_phy_sf_vars_setup(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver(), and mod_atmos_driver::atmos_surface_get().

67  real(RP), public, allocatable :: atmos_phy_sf_sflx_mu (:,:) ! x-momentum flux (area center) [m/s*kg/m2/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu

◆ atmos_phy_sf_sflx_mv

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv

Definition at line 68 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup(), atmos_phy_sf_vars_setup(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver(), and mod_atmos_driver::atmos_surface_get().

68  real(RP), public, allocatable :: atmos_phy_sf_sflx_mv (:,:) ! y-momentum flux (area center) [m/s*kg/m2/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv

◆ atmos_phy_sf_sflx_sh

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh

◆ atmos_phy_sf_sflx_lh

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_lh

Definition at line 70 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup(), atmos_phy_sf_vars_setup(), mod_atmos_driver::atmos_surface_get(), and mod_atmos_vars::atmos_vars_monitor().

70  real(RP), public, allocatable :: atmos_phy_sf_sflx_lh (:,:) ! latent heat flux [J/m2/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh

◆ atmos_phy_sf_sflx_gh

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_gh

Definition at line 71 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), atmos_phy_sf_vars_setup(), and mod_atmos_driver::atmos_surface_get().

71  real(RP), public, allocatable :: atmos_phy_sf_sflx_gh (:,:) ! ground heat flux [J/m2/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh

◆ atmos_phy_sf_sflx_qtrc

real(rp), dimension (:,:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc

◆ atmos_phy_sf_u10

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_u10

Definition at line 74 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), atmos_phy_sf_vars_setup(), and mod_atmos_driver::atmos_surface_get().

74  real(RP), public, allocatable :: atmos_phy_sf_u10 (:,:) ! 10m x-wind [m/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10

◆ atmos_phy_sf_v10

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_v10

Definition at line 75 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), atmos_phy_sf_vars_setup(), and mod_atmos_driver::atmos_surface_get().

75  real(RP), public, allocatable :: atmos_phy_sf_v10 (:,:) ! 10m y-wind [m/s]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_v10

◆ atmos_phy_sf_t2

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_t2

Definition at line 76 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), atmos_phy_sf_vars_setup(), and mod_atmos_driver::atmos_surface_get().

76  real(RP), public, allocatable :: atmos_phy_sf_t2 (:,:) ! 2m temperature [K]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2

◆ atmos_phy_sf_q2

real(rp), dimension (:,:), allocatable, public mod_atmos_phy_sf_vars::atmos_phy_sf_q2

Definition at line 77 of file mod_atmos_phy_sf_vars.f90.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), atmos_phy_sf_vars_setup(), and mod_atmos_driver::atmos_surface_get().

77  real(RP), public, allocatable :: atmos_phy_sf_q2 (:,:) ! 2m specific humidity [kg/kg]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2