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_open
 Open restart file for read. More...
 
subroutine, public atmos_phy_sf_vars_restart_read
 Read 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...
 
subroutine, public atmos_phy_sf_vars_restart_create
 Create restart file. More...
 
subroutine, public atmos_phy_sf_vars_restart_enddef
 Exit netCDF define mode. More...
 
subroutine, public atmos_phy_sf_vars_restart_close
 Close restart file. More...
 
subroutine, public atmos_phy_sf_vars_restart_def_var
 Write restart. More...
 
subroutine, public atmos_phy_sf_vars_restart_write
 Write variables to restart file. 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 input file. More...
 
logical, public atmos_phy_sf_restart_in_postfix_timelabel = .false.
 Add timelabel to the basename of input file? More...
 
character(len=h_long), public atmos_phy_sf_restart_out_basename = ''
 Basename of the output file. More...
 
logical, public atmos_phy_sf_restart_out_postfix_timelabel = .true.
 Add timelabel to the basename of 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_short), 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 input file
    ATMOS_PHY_SF_RESTART_IN_POSTFIX_TIMELABEL logical .false. Add timelabel to the basename of input 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_POSTFIX_TIMELABEL logical .true. Add timelabel to the basename of 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_SHORT) '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 143 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_in_postfix_timelabel, atmos_phy_sf_restart_out_basename, atmos_phy_sf_restart_out_dtype, atmos_phy_sf_restart_out_postfix_timelabel, 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_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_grid_index::ja, scale_process::prc_mpistop(), and scale_tracer::qa.

Referenced by mod_atmos_vars::atmos_vars_setup().

143  use scale_process, only: &
145  use scale_const, only: &
146  undef => const_undef
147  implicit none
148 
149  namelist / param_atmos_phy_sf_vars / &
150  atmos_phy_sf_restart_in_basename, &
151  atmos_phy_sf_restart_in_postfix_timelabel, &
152  atmos_phy_sf_restart_output, &
153  atmos_phy_sf_restart_out_basename, &
154  atmos_phy_sf_restart_out_postfix_timelabel, &
155  atmos_phy_sf_restart_out_title, &
156  atmos_phy_sf_restart_out_dtype, &
157  atmos_phy_sf_default_sfc_temp, &
158  atmos_phy_sf_default_sfc_albedo, &
159  atmos_phy_sf_default_sfc_z0
160 
161  integer :: ierr
162  integer :: iv
163  !---------------------------------------------------------------------------
164 
165  if( io_l ) write(io_fid_log,*)
166  if( io_l ) write(io_fid_log,*) '++++++ Module[VARS] / Categ[ATMOS PHY_SF] / Origin[SCALE-RM]'
167 
168  allocate( atmos_phy_sf_dens_t(ia,ja) )
169  allocate( atmos_phy_sf_momz_t(ia,ja) )
170  allocate( atmos_phy_sf_momx_t(ia,ja) )
171  allocate( atmos_phy_sf_momy_t(ia,ja) )
172  allocate( atmos_phy_sf_rhot_t(ia,ja) )
173  allocate( atmos_phy_sf_rhoq_t(ia,ja,qa) )
174  atmos_phy_sf_dens_t(:,:) = undef
175  atmos_phy_sf_momz_t(:,:) = undef
176  atmos_phy_sf_momx_t(:,:) = undef
177  atmos_phy_sf_momy_t(:,:) = undef
178  atmos_phy_sf_rhot_t(:,:) = undef
179  atmos_phy_sf_rhoq_t(:,:,:) = undef
180 
181  allocate( atmos_phy_sf_sfc_temp(ia,ja) )
182  allocate( atmos_phy_sf_sfc_albedo(ia,ja,2) )
183  allocate( atmos_phy_sf_sfc_z0m(ia,ja) )
184  allocate( atmos_phy_sf_sfc_z0h(ia,ja) )
185  allocate( atmos_phy_sf_sfc_z0e(ia,ja) )
186 ! ATMOS_PHY_SF_SFC_TEMP (:,:) = UNDEF ! [del] 2014/8/28 A.Noda
187 ! ATMOS_PHY_SF_SFC_albedo(:,:,:) = UNDEF ! [del] 2014/8/28 A.Noda
188  atmos_phy_sf_sfc_z0m(:,:) = undef
189  atmos_phy_sf_sfc_z0h(:,:) = undef
190  atmos_phy_sf_sfc_z0e(:,:) = undef
191 
192  allocate( atmos_phy_sf_sfc_dens(ia,ja) )
193  allocate( atmos_phy_sf_sfc_pres(ia,ja) )
194  atmos_phy_sf_sfc_dens(:,:) = undef
195  atmos_phy_sf_sfc_pres(:,:) = undef
196 
197  allocate( atmos_phy_sf_sflx_mw(ia,ja) )
198  allocate( atmos_phy_sf_sflx_mu(ia,ja) )
199  allocate( atmos_phy_sf_sflx_mv(ia,ja) )
200  allocate( atmos_phy_sf_sflx_sh(ia,ja) )
201  allocate( atmos_phy_sf_sflx_lh(ia,ja) )
202  allocate( atmos_phy_sf_sflx_gh(ia,ja) )
203  allocate( atmos_phy_sf_sflx_qtrc(ia,ja,qa) )
204  atmos_phy_sf_sflx_mw(:,:) = undef
205  atmos_phy_sf_sflx_mu(:,:) = undef
206  atmos_phy_sf_sflx_mv(:,:) = undef
207  atmos_phy_sf_sflx_sh(:,:) = undef
208  atmos_phy_sf_sflx_lh(:,:) = undef
209  atmos_phy_sf_sflx_gh(:,:) = undef
210  atmos_phy_sf_sflx_qtrc(:,:,:) = undef
211 
212  allocate( atmos_phy_sf_u10(ia,ja) )
213  allocate( atmos_phy_sf_v10(ia,ja) )
214  allocate( atmos_phy_sf_t2(ia,ja) )
215  allocate( atmos_phy_sf_q2(ia,ja) )
216  atmos_phy_sf_u10(:,:) = undef
217  atmos_phy_sf_v10(:,:) = undef
218  atmos_phy_sf_t2(:,:) = undef
219  atmos_phy_sf_q2(:,:) = undef
220 
221  !--- read namelist
222  rewind(io_fid_conf)
223  read(io_fid_conf,nml=param_atmos_phy_sf_vars,iostat=ierr)
224  if( ierr < 0 ) then !--- missing
225  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
226  elseif( ierr > 0 ) then !--- fatal error
227  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_VARS. Check!'
228  call prc_mpistop
229  endif
230  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_sf_vars)
231 
232  ! [add] 2014/08/28 A.Noda
233  atmos_phy_sf_sfc_temp(:,:) = atmos_phy_sf_default_sfc_temp
234  atmos_phy_sf_sfc_albedo(:,:,:) = atmos_phy_sf_default_sfc_albedo
235  atmos_phy_sf_sfc_z0m(:,:) = atmos_phy_sf_default_sfc_z0
236  atmos_phy_sf_sfc_z0h(:,:) = atmos_phy_sf_default_sfc_z0
237  atmos_phy_sf_sfc_z0e(:,:) = atmos_phy_sf_default_sfc_z0
238 
239  if( io_l ) write(io_fid_log,*)
240  if( io_l ) write(io_fid_log,*) '*** [ATMOS_PHY_SF] prognostic/diagnostic variables'
241  if( io_l ) write(io_fid_log,'(1x,A,A24,A,A48,A,A12,A)') &
242  '*** |', 'VARNAME ','|', &
243  'DESCRIPTION ', '[', 'UNIT ', ']'
244  do iv = 1, vmax
245  if( io_l ) write(io_fid_log,'(1x,A,I3,A,A24,A,A48,A,A12,A)') &
246  '*** NO.',iv,'|',var_name(iv),'|',var_desc(iv),'[',var_unit(iv),']'
247  enddo
248 
249  if( io_l ) write(io_fid_log,*)
250  if ( atmos_phy_sf_restart_in_basename /= '' ) then
251  if( io_l ) write(io_fid_log,*) '*** Restart input? : YES, file = ', trim(atmos_phy_sf_restart_in_basename)
252  if( io_l ) write(io_fid_log,*) '*** Add timelabel? : ', atmos_phy_sf_restart_in_postfix_timelabel
253  else
254  if( io_l ) write(io_fid_log,*) '*** Restart input? : NO'
255  endif
256  if ( atmos_phy_sf_restart_output &
257  .AND. atmos_phy_sf_restart_out_basename /= '' ) then
258  if( io_l ) write(io_fid_log,*) '*** Restart output? : YES, file = ', trim(atmos_phy_sf_restart_out_basename)
259  if( io_l ) write(io_fid_log,*) '*** Add timelabel? : ', atmos_phy_sf_restart_out_postfix_timelabel
260  else
261  if( io_l ) write(io_fid_log,*) '*** Restart output? : NO'
262  atmos_phy_sf_restart_output = .false.
263  endif
264 
265  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
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhot_t
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
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
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_fillhalo()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_fillhalo ( )

HALO Communication.

Definition at line 271 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().

271  use scale_const, only: &
272  i_sw => const_i_sw, &
273  i_lw => const_i_lw
274  use scale_comm, only: &
275  comm_vars8, &
276  comm_wait
277  implicit none
278  !---------------------------------------------------------------------------
279 
280  call comm_vars8( atmos_phy_sf_sfc_temp(:,:), 1 )
281  call comm_vars8( atmos_phy_sf_sfc_albedo(:,:,i_lw), 2 )
282  call comm_vars8( atmos_phy_sf_sfc_albedo(:,:,i_sw), 3 )
283  call comm_vars8( atmos_phy_sf_sfc_z0m(:,:), 4 )
284  call comm_vars8( atmos_phy_sf_sfc_z0h(:,:), 5 )
285  call comm_vars8( atmos_phy_sf_sfc_z0e(:,:), 6 )
286 
287  call comm_wait ( atmos_phy_sf_sfc_temp(:,:), 1 )
288  call comm_wait ( atmos_phy_sf_sfc_albedo(:,:,i_lw), 2 )
289  call comm_wait ( atmos_phy_sf_sfc_albedo(:,:,i_sw), 3 )
290  call comm_wait ( atmos_phy_sf_sfc_z0m(:,:), 4 )
291  call comm_wait ( atmos_phy_sf_sfc_z0h(:,:), 5 )
292  call comm_wait ( atmos_phy_sf_sfc_z0e(:,:), 6 )
293 
294  return
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:95
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:96
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_restart_open()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_open ( )

Open restart file for read.

Definition at line 300 of file mod_atmos_phy_sf_vars.f90.

References atmos_phy_sf_restart_in_basename, atmos_phy_sf_restart_in_postfix_timelabel, 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_fileio::fileio_open(), scale_stdio::io_fid_log, scale_stdio::io_l, and scale_time::time_gettimelabel().

Referenced by mod_atmos_vars::atmos_vars_restart_open().

300  use scale_time, only: &
302  use scale_fileio, only: &
304  implicit none
305 
306  character(len=19) :: timelabel
307  character(len=H_LONG) :: basename
308  !---------------------------------------------------------------------------
309 
310  if( io_l ) write(io_fid_log,*)
311  if( io_l ) write(io_fid_log,*) '*** Open restart file (ATMOS_PHY_SF) ***'
312 
313  if ( atmos_phy_sf_restart_in_basename /= '' ) then
314 
315  if ( atmos_phy_sf_restart_in_postfix_timelabel ) then
316  call time_gettimelabel( timelabel )
317  basename = trim(atmos_phy_sf_restart_in_basename)//'_'//trim(timelabel)
318  else
319  basename = trim(atmos_phy_sf_restart_in_basename)
320  endif
321 
322  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(basename)
323 
324  call fileio_open( restart_fid, basename )
325 
326  else
327  if( io_l ) write(io_fid_log,*) '*** restart file for ATMOS_PHY_SF is not specified.'
328  atmos_phy_sf_sfc_temp(:,:) = atmos_phy_sf_default_sfc_temp
329  atmos_phy_sf_sfc_albedo(:,:,:) = atmos_phy_sf_default_sfc_albedo
330  atmos_phy_sf_sfc_z0m(:,:) = atmos_phy_sf_default_sfc_z0
331  atmos_phy_sf_sfc_z0h(:,:) = atmos_phy_sf_default_sfc_z0
332  atmos_phy_sf_sfc_z0e(:,:) = atmos_phy_sf_default_sfc_z0
333  endif
334 
335  return
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
module FILE I/O (netcdf)
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
subroutine, public fileio_open(fid, basename)
open a netCDF file for read
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
Here is the call graph for this function:
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 341 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_fileio::fileio_flush(), scale_stdio::io_aggregate, scale_stdio::io_fid_log, scale_stdio::io_l, and scale_rm_statistics::statistics_checktotal.

Referenced by mod_atmos_vars::atmos_vars_restart_read().

341  use scale_const, only: &
342  i_sw => const_i_sw, &
343  i_lw => const_i_lw
344  use scale_rm_statistics, only: &
346  stat_total
347  use scale_fileio, only: &
348  fileio_read, &
350  implicit none
351 
352  real(RP) :: total
353  !---------------------------------------------------------------------------
354 
355  if ( restart_fid /= -1 ) then
356  if( io_l ) write(io_fid_log,*)
357  if( io_l ) write(io_fid_log,*) '*** Read from restart file (ATMOS_PHY_SF) ***'
358 
359  call fileio_read( atmos_phy_sf_sfc_temp(:,:), & ! [OUT]
360  restart_fid, var_name(1), 'XY', step=1 ) ! [IN]
361  call fileio_read( atmos_phy_sf_sfc_albedo(:,:,i_lw), & ! [OUT]
362  restart_fid, var_name(2), 'XY', step=1 ) ! [IN]
363  call fileio_read( atmos_phy_sf_sfc_albedo(:,:,i_sw), & ! [OUT]
364  restart_fid, var_name(3), 'XY', step=1 ) ! [IN]
365  call fileio_read( atmos_phy_sf_sfc_z0m(:,:), & ! [OUT]
366  restart_fid, var_name(4), 'XY', step=1 ) ! [IN]
367  call fileio_read( atmos_phy_sf_sfc_z0h(:,:), & ! [OUT]
368  restart_fid, var_name(5), 'XY', step=1 ) ! [IN]
369  call fileio_read( atmos_phy_sf_sfc_z0e(:,:), & ! [OUT]
370  restart_fid, var_name(6), 'XY', step=1 ) ! [IN]
371 
372  if ( io_aggregate) then
373  call fileio_flush( restart_fid ) ! X/Y halos have been read from file
374  else
375  call atmos_phy_sf_vars_fillhalo
376  end if
377 
378  if ( statistics_checktotal ) then
379  call stat_total( total, atmos_phy_sf_sfc_temp(:,:), var_name(1) )
380  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_lw), var_name(2) )
381  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_sw), var_name(3) )
382  call stat_total( total, atmos_phy_sf_sfc_z0m(:,:), var_name(4) )
383  call stat_total( total, atmos_phy_sf_sfc_z0h(:,:), var_name(5) )
384  call stat_total( total, atmos_phy_sf_sfc_z0e(:,:), var_name(6) )
385  end if
386  else
387  if( io_l ) write(io_fid_log,*) '*** invalid restart file ID for ATMOS_PHY_SF.'
388  endif
389 
390  return
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:95
subroutine, public fileio_flush(fid)
Flush all pending requests to a netCDF file (PnetCDF only)
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:96
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
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 402 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.

402  use scale_const, only: &
403  i_sw => const_i_sw, &
404  i_lw => const_i_lw
405  implicit none
406 
407  real(RP), intent(in) :: sfc_temp_in (:,:)
408  real(RP), intent(in) :: albedo_lw_in(:,:)
409  real(RP), intent(in) :: albedo_sw_in(:,:)
410  real(RP), intent(in) :: sfc_z0m_in (:,:)
411  real(RP), intent(in) :: sfc_z0h_in (:,:)
412  real(RP), intent(in) :: sfc_z0e_in (:,:)
413  !---------------------------------------------------------------------------
414 
415  if( io_l ) write(io_fid_log,*)
416  if( io_l ) write(io_fid_log,*) '*** External Input (PHY_SF) ***'
417 
418  atmos_phy_sf_sfc_temp(:,:) = sfc_temp_in(:,:)
419  atmos_phy_sf_sfc_albedo(:,:,i_lw) = albedo_lw_in(:,:)
420  atmos_phy_sf_sfc_albedo(:,:,i_sw) = albedo_sw_in(:,:)
421  atmos_phy_sf_sfc_z0m(:,:) = sfc_z0m_in(:,:)
422  atmos_phy_sf_sfc_z0h(:,:) = sfc_z0h_in(:,:)
423  atmos_phy_sf_sfc_z0e(:,:) = sfc_z0e_in(:,:)
424 
425  call atmos_phy_sf_vars_fillhalo
426 
427  return
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:95
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:96
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
Here is the call graph for this function:

◆ atmos_phy_sf_vars_restart_create()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_create ( )

Create restart file.

Definition at line 433 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_postfix_timelabel, atmos_phy_sf_restart_out_title, scale_fileio::fileio_create(), scale_stdio::io_fid_log, scale_stdio::io_l, and scale_time::time_gettimelabel().

Referenced by mod_atmos_vars::atmos_vars_restart_create().

433  use scale_time, only: &
435  use scale_fileio, only: &
437  implicit none
438 
439  character(len=19) :: timelabel
440  character(len=H_LONG) :: basename
441  !---------------------------------------------------------------------------
442 
443  if ( atmos_phy_sf_restart_out_basename /= '' ) then
444 
445  if( io_l ) write(io_fid_log,*)
446  if( io_l ) write(io_fid_log,*) '*** Create restart file (ATMOS_PHY_AE) ***'
447 
448  if ( atmos_phy_sf_restart_out_postfix_timelabel ) then
449  call time_gettimelabel( timelabel )
450  basename = trim(atmos_phy_sf_restart_out_basename)//'_'//trim(timelabel)
451  else
452  basename = trim(atmos_phy_sf_restart_out_basename)
453  endif
454 
455  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(basename)
456 
457  call fileio_create( restart_fid, & ! [OUT]
458  basename, atmos_phy_sf_restart_out_title, atmos_phy_sf_restart_out_dtype ) ! [IN]
459  endif
460 
461  return
module FILE I/O (netcdf)
subroutine, public time_gettimelabel(timelabel)
generate time label
Definition: scale_time.F90:90
subroutine, public fileio_create(fid, basename, title, datatype, date, subsec, append, nozcoord)
Create/open a netCDF file.
module TIME
Definition: scale_time.F90:15
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_restart_enddef()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_enddef ( )

Exit netCDF define mode.

Definition at line 467 of file mod_atmos_phy_sf_vars.f90.

References scale_fileio::fileio_enddef().

Referenced by mod_atmos_vars::atmos_vars_restart_enddef().

467  use scale_fileio, only: &
469  implicit none
470 
471  if ( restart_fid /= -1 ) then
472  call fileio_enddef( restart_fid ) ! [IN]
473  endif
474 
475  return
module FILE I/O (netcdf)
subroutine, public fileio_enddef(fid)
Exit netCDF file define mode.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_restart_close()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_close ( )

Close restart file.

Definition at line 481 of file mod_atmos_phy_sf_vars.f90.

References scale_fileio::fileio_close(), scale_stdio::io_fid_log, and scale_stdio::io_l.

Referenced by mod_atmos_vars::atmos_vars_restart_close().

481  use scale_fileio, only: &
483  implicit none
484  !---------------------------------------------------------------------------
485 
486  if ( restart_fid /= -1 ) then
487  if( io_l ) write(io_fid_log,*)
488  if( io_l ) write(io_fid_log,*) '*** Close restart file (ATMOS_PHY_SF) ***'
489 
490  call fileio_close( restart_fid ) ! [IN]
491 
492  restart_fid = -1
493  endif
494 
495  return
module FILE I/O (netcdf)
subroutine, public fileio_close(fid)
Close a netCDF file.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_vars_restart_def_var()

subroutine, public mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_def_var ( )

Write restart.

Definition at line 501 of file mod_atmos_phy_sf_vars.f90.

References atmos_phy_sf_restart_out_dtype, and scale_fileio::fileio_def_var().

Referenced by mod_atmos_vars::atmos_vars_restart_def_var().

501  use scale_fileio, only: &
503  implicit none
504  !---------------------------------------------------------------------------
505 
506  if ( restart_fid /= -1 ) then
507 
508  call fileio_def_var( restart_fid, var_id(1), var_name(1), var_desc(1), &
509  var_unit(1), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
510  call fileio_def_var( restart_fid, var_id(2), var_name(2), var_desc(2), &
511  var_unit(2), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
512  call fileio_def_var( restart_fid, var_id(3), var_name(3), var_desc(3), &
513  var_unit(3), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
514  call fileio_def_var( restart_fid, var_id(4), var_name(4), var_desc(4), &
515  var_unit(4), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
516  call fileio_def_var( restart_fid, var_id(5), var_name(5), var_desc(5), &
517  var_unit(5), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
518  call fileio_def_var( restart_fid, var_id(6), var_name(6), var_desc(6), &
519  var_unit(6), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
520  endif
521 
522  return
module FILE I/O (netcdf)
subroutine, public fileio_def_var(fid, vid, varname, desc, unit, axistype, datatype, timeintv, nsteps)
Define a variable to file.
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 variables to restart file.

Definition at line 528 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, and scale_rm_statistics::statistics_checktotal.

Referenced by mod_atmos_vars::atmos_vars_restart_write().

528  use scale_const, only: &
529  i_sw => const_i_sw, &
530  i_lw => const_i_lw
531  use scale_rm_statistics, only: &
533  stat_total
534  use scale_fileio, only: &
535  fileio_write_var
536  implicit none
537 
538  real(RP) :: total
539  !---------------------------------------------------------------------------
540 
541  if ( restart_fid /= -1 ) then
542 
543  call atmos_phy_sf_vars_fillhalo
544 
545  if ( statistics_checktotal ) then
546  call stat_total( total, atmos_phy_sf_sfc_temp(:,:), var_name(1) )
547  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_lw), var_name(2) )
548  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_sw), var_name(3) )
549  call stat_total( total, atmos_phy_sf_sfc_z0m(:,:), var_name(4) )
550  call stat_total( total, atmos_phy_sf_sfc_z0h(:,:), var_name(5) )
551  call stat_total( total, atmos_phy_sf_sfc_z0e(:,:), var_name(6) )
552  endif
553 
554  call fileio_write_var( restart_fid, var_id(1), atmos_phy_sf_sfc_temp(:,:), & ! [IN]
555  var_name(1), 'XY' ) ! [IN]
556  call fileio_write_var( restart_fid, var_id(2), atmos_phy_sf_sfc_albedo(:,:,i_lw), & ! [IN]
557  var_name(2), 'XY' ) ! [IN]
558  call fileio_write_var( restart_fid, var_id(3), atmos_phy_sf_sfc_albedo(:,:,i_sw), & ! [IN]
559  var_name(3), 'XY' ) ! [IN]
560  call fileio_write_var( restart_fid, var_id(4), atmos_phy_sf_sfc_z0m(:,:), & ! [IN]
561  var_name(4), 'XY' ) ! [IN]
562  call fileio_write_var( restart_fid, var_id(5), atmos_phy_sf_sfc_z0h(:,:), & ! [IN]
563  var_name(5), 'XY' ) ! [IN]
564  call fileio_write_var( restart_fid, var_id(6), atmos_phy_sf_sfc_z0e(:,:), & ! [IN]
565  var_name(6), 'XY' ) ! [IN]
566 
567  endif
568 
569  return
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:95
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:96
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
Here is the call graph for this function:
Here is the caller 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 48 of file mod_atmos_phy_sf_vars.f90.

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

48  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 input file.

Definition at line 50 of file mod_atmos_phy_sf_vars.f90.

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

50  character(len=H_LONG), public :: ATMOS_PHY_SF_RESTART_IN_BASENAME = ''

◆ atmos_phy_sf_restart_in_postfix_timelabel

logical, public mod_atmos_phy_sf_vars::atmos_phy_sf_restart_in_postfix_timelabel = .false.

Add timelabel to the basename of input file?

Definition at line 51 of file mod_atmos_phy_sf_vars.f90.

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

51  logical, public :: ATMOS_PHY_SF_RESTART_IN_POSTFIX_TIMELABEL = .false.

◆ 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 52 of file mod_atmos_phy_sf_vars.f90.

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

52  character(len=H_LONG), public :: ATMOS_PHY_SF_RESTART_OUT_BASENAME = ''

◆ atmos_phy_sf_restart_out_postfix_timelabel

logical, public mod_atmos_phy_sf_vars::atmos_phy_sf_restart_out_postfix_timelabel = .true.

Add timelabel to the basename of output file?

Definition at line 53 of file mod_atmos_phy_sf_vars.f90.

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

53  logical, public :: ATMOS_PHY_SF_RESTART_OUT_POSTFIX_TIMELABEL = .true.

◆ 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 54 of file mod_atmos_phy_sf_vars.f90.

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

54  character(len=H_MID), public :: ATMOS_PHY_SF_RESTART_OUT_TITLE = 'ATMOS_PHY_SF restart'

◆ atmos_phy_sf_restart_out_dtype

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

REAL4 or REAL8.

Definition at line 55 of file mod_atmos_phy_sf_vars.f90.

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

55  character(len=H_SHORT), 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 57 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().

57  real(RP), public, allocatable :: ATMOS_PHY_SF_DENS_t(:,:) ! tendency DENS [ kg/m3/s]

◆ atmos_phy_sf_momz_t

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

Definition at line 58 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().

58  real(RP), public, allocatable :: ATMOS_PHY_SF_MOMZ_t(:,:) ! tendency MOMZ [m/s*kg/m3/s]

◆ atmos_phy_sf_momx_t

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

Definition at line 59 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().

59  real(RP), public, allocatable :: ATMOS_PHY_SF_MOMX_t(:,:) ! tendency MOMX [m/s*kg/m3/s]

◆ atmos_phy_sf_momy_t

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

Definition at line 60 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().

60  real(RP), public, allocatable :: ATMOS_PHY_SF_MOMY_t(:,:) ! tendency MOMY [m/s*kg/m3/s]

◆ atmos_phy_sf_rhot_t

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

Definition at line 61 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().

61  real(RP), public, allocatable :: ATMOS_PHY_SF_RHOT_t(:,:) ! tendency RHOT [K *kg/m3/s]

◆ atmos_phy_sf_rhoq_t

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

Definition at line 62 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().

62  real(RP), public, allocatable :: ATMOS_PHY_SF_RHOQ_t(:,:,:) ! tendency rho*QTRC [ kg/kg/s]

◆ 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 71 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().

71  real(RP), public, allocatable :: ATMOS_PHY_SF_SFC_DENS (:,:) ! surface atmosphere density [kg/m3]

◆ atmos_phy_sf_sfc_pres

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

Definition at line 72 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().

72  real(RP), public, allocatable :: ATMOS_PHY_SF_SFC_PRES (:,:) ! surface atmosphere pressure [Pa]

◆ atmos_phy_sf_sflx_mw

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

Definition at line 74 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().

74  real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_MW (:,:) ! z-momentum flux (area center) [m/s*kg/m2/s]

◆ atmos_phy_sf_sflx_mu

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

Definition at line 75 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().

75  real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_MU (:,:) ! x-momentum flux (area center) [m/s*kg/m2/s]

◆ atmos_phy_sf_sflx_mv

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

Definition at line 76 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().

76  real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_MV (:,:) ! y-momentum flux (area center) [m/s*kg/m2/s]

◆ 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 78 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().

78  real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_LH (:,:) ! latent heat flux [J/m2/s]

◆ atmos_phy_sf_sflx_gh

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

Definition at line 79 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().

79  real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_GH (:,:) ! ground heat flux [J/m2/s]

◆ 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 82 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().

82  real(RP), public, allocatable :: ATMOS_PHY_SF_U10 (:,:) ! 10m x-wind [m/s]

◆ atmos_phy_sf_v10

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

Definition at line 83 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().

83  real(RP), public, allocatable :: ATMOS_PHY_SF_V10 (:,:) ! 10m y-wind [m/s]

◆ atmos_phy_sf_t2

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

Definition at line 84 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().

84  real(RP), public, allocatable :: ATMOS_PHY_SF_T2 (:,:) ! 2m temperature [K]

◆ atmos_phy_sf_q2

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

Definition at line 85 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().

85  real(RP), public, allocatable :: ATMOS_PHY_SF_Q2 (:,:) ! 2m specific humidity [kg/kg]