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 
43 
44  !-----------------------------------------------------------------------------
45  !
46  !++ Public parameters & variables
47  !
48  logical, public :: atmos_phy_sf_restart_output = .false.
49 
50  character(len=H_LONG), public :: atmos_phy_sf_restart_in_basename = ''
51  character(len=H_LONG), public :: atmos_phy_sf_restart_out_basename = ''
52  character(len=H_MID), public :: atmos_phy_sf_restart_out_title = 'ATMOS_PHY_SF restart'
53  character(len=H_MID), public :: atmos_phy_sf_restart_out_dtype = 'DEFAULT'
54 
55  real(RP), public, allocatable :: atmos_phy_sf_dens_t(:,:) ! tendency DENS [ kg/m3/s]
56  real(RP), public, allocatable :: atmos_phy_sf_momz_t(:,:) ! tendency MOMZ [m/s*kg/m3/s]
57  real(RP), public, allocatable :: atmos_phy_sf_momx_t(:,:) ! tendency MOMX [m/s*kg/m3/s]
58  real(RP), public, allocatable :: atmos_phy_sf_momy_t(:,:) ! tendency MOMY [m/s*kg/m3/s]
59  real(RP), public, allocatable :: atmos_phy_sf_rhot_t(:,:) ! tendency RHOT [K *kg/m3/s]
60  real(RP), public, allocatable :: atmos_phy_sf_rhoq_t(:,:,:) ! tendency rho*QTRC [ kg/kg/s]
61 
62  real(RP), public, allocatable :: atmos_phy_sf_sfc_temp (:,:) ! surface skin temperature [K]
63  real(RP), public, allocatable :: atmos_phy_sf_sfc_albedo(:,:,:) ! surface albedo [0-1]
64  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0m (:,:) ! surface roughness length, ocean only [m]
65  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0h (:,:) ! surface roughness length, ocean only [m]
66  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0e (:,:) ! surface roughness length, ocean only [m]
67 
68  ! surface diagnostic variables
69  real(RP), public, allocatable :: atmos_phy_sf_sfc_dens (:,:) ! surface atmosphere density [kg/m3]
70  real(RP), public, allocatable :: atmos_phy_sf_sfc_pres (:,:) ! surface atmosphere pressure [Pa]
71 
72  real(RP), public, allocatable :: atmos_phy_sf_sflx_mw (:,:) ! z-momentum flux (area center) [m/s*kg/m2/s]
73  real(RP), public, allocatable :: atmos_phy_sf_sflx_mu (:,:) ! x-momentum flux (area center) [m/s*kg/m2/s]
74  real(RP), public, allocatable :: atmos_phy_sf_sflx_mv (:,:) ! y-momentum flux (area center) [m/s*kg/m2/s]
75  real(RP), public, allocatable :: atmos_phy_sf_sflx_sh (:,:) ! sensible heat flux [J/m2/s]
76  real(RP), public, allocatable :: atmos_phy_sf_sflx_lh (:,:) ! latent heat flux [J/m2/s]
77  real(RP), public, allocatable :: atmos_phy_sf_sflx_gh (:,:) ! ground heat flux [J/m2/s]
78  real(RP), public, allocatable :: atmos_phy_sf_sflx_qtrc (:,:,:) ! tracer mass flux [kg/m2/s]
79 
80  real(RP), public, allocatable :: atmos_phy_sf_u10 (:,:) ! 10m x-wind [m/s]
81  real(RP), public, allocatable :: atmos_phy_sf_v10 (:,:) ! 10m y-wind [m/s]
82  real(RP), public, allocatable :: atmos_phy_sf_t2 (:,:) ! 2m temperature [K]
83  real(RP), public, allocatable :: atmos_phy_sf_q2 (:,:) ! 2m specific humidity [kg/kg]
84 
85 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_QEMIS(:,:,:) ! tracer emission flux [kg/m2/s]
86 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_QDEP (:,:,:) ! tracer deposition flux [kg/m2/s]
87 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_VDEP (:,:,:) ! tracer deposition velocity [m/s]
88 
89  !-----------------------------------------------------------------------------
90  !
91  !++ Private procedure
92  !
93  !-----------------------------------------------------------------------------
94  !
95  !++ Private parameters & variables
96  !
97  integer, private, parameter :: vmax = 6
98  integer, private, parameter :: i_sfc_temp = 1
99  integer, private, parameter :: i_sfc_albedo_lw = 2
100  integer, private, parameter :: i_sfc_albedo_sw = 3
101  integer, private, parameter :: i_sfc_z0m = 4
102  integer, private, parameter :: i_sfc_z0h = 5
103  integer, private, parameter :: i_sfc_z0e = 6
104 
105  character(len=H_SHORT), private :: var_name(vmax)
106  character(len=H_MID), private :: var_desc(vmax)
107  character(len=H_SHORT), private :: var_unit(vmax)
108  integer, private :: var_id(vmax)
109  integer, private :: restart_fid = -1 ! file ID
110 
111  data var_name / 'SFC_TEMP', &
112  'SFC_ALB_LW', &
113  'SFC_ALB_SW', &
114  'SFC_Z0M', &
115  'SFC_Z0H', &
116  'SFC_Z0E' /
117 
118  data var_desc / 'surface skin temperature', &
119  'surface albedo for longwave', &
120  'surface albedo for shortwave', &
121  'surface roughness length (momentum)', &
122  'surface roughness length (heat)', &
123  'surface roughness length (vapor)' /
124 
125  data var_unit / 'K', &
126  '0-1', &
127  '0-1', &
128  'm', &
129  'm', &
130  'm' /
131 
132  real(RP), private :: atmos_phy_sf_default_sfc_temp = 300.0_rp
133  real(RP), private :: atmos_phy_sf_default_sfc_albedo = 0.4_rp
134  real(RP), private :: atmos_phy_sf_default_sfc_z0 = 1e-4_rp
135 
136  !-----------------------------------------------------------------------------
137 contains
138  !-----------------------------------------------------------------------------
140  subroutine atmos_phy_sf_vars_setup
141  use scale_process, only: &
143  use scale_const, only: &
144  undef => const_undef
145  implicit none
146 
147  namelist / param_atmos_phy_sf_vars / &
153  atmos_phy_sf_default_sfc_temp, &
154  atmos_phy_sf_default_sfc_albedo, &
155  atmos_phy_sf_default_sfc_z0
156 
157  integer :: ierr
158  integer :: iv
159  !---------------------------------------------------------------------------
160 
161  if( io_l ) write(io_fid_log,*)
162  if( io_l ) write(io_fid_log,*) '++++++ Module[VARS] / Categ[ATMOS PHY_SF] / Origin[SCALE-RM]'
163 
164  allocate( atmos_phy_sf_dens_t(ia,ja) )
165  allocate( atmos_phy_sf_momz_t(ia,ja) )
166  allocate( atmos_phy_sf_momx_t(ia,ja) )
167  allocate( atmos_phy_sf_momy_t(ia,ja) )
168  allocate( atmos_phy_sf_rhot_t(ia,ja) )
169  allocate( atmos_phy_sf_rhoq_t(ia,ja,qa) )
170  atmos_phy_sf_dens_t(:,:) = undef
171  atmos_phy_sf_momz_t(:,:) = undef
172  atmos_phy_sf_momx_t(:,:) = undef
173  atmos_phy_sf_momy_t(:,:) = undef
174  atmos_phy_sf_rhot_t(:,:) = undef
175  atmos_phy_sf_rhoq_t(:,:,:) = undef
176 
177  allocate( atmos_phy_sf_sfc_temp(ia,ja) )
178  allocate( atmos_phy_sf_sfc_albedo(ia,ja,2) )
179  allocate( atmos_phy_sf_sfc_z0m(ia,ja) )
180  allocate( atmos_phy_sf_sfc_z0h(ia,ja) )
181  allocate( atmos_phy_sf_sfc_z0e(ia,ja) )
182 ! ATMOS_PHY_SF_SFC_TEMP (:,:) = UNDEF ! [del] 2014/8/28 A.Noda
183 ! ATMOS_PHY_SF_SFC_albedo(:,:,:) = UNDEF ! [del] 2014/8/28 A.Noda
184  atmos_phy_sf_sfc_z0m(:,:) = undef
185  atmos_phy_sf_sfc_z0h(:,:) = undef
186  atmos_phy_sf_sfc_z0e(:,:) = undef
187 
188  allocate( atmos_phy_sf_sfc_dens(ia,ja) )
189  allocate( atmos_phy_sf_sfc_pres(ia,ja) )
190  atmos_phy_sf_sfc_dens(:,:) = undef
191  atmos_phy_sf_sfc_pres(:,:) = undef
192 
193  allocate( atmos_phy_sf_sflx_mw(ia,ja) )
194  allocate( atmos_phy_sf_sflx_mu(ia,ja) )
195  allocate( atmos_phy_sf_sflx_mv(ia,ja) )
196  allocate( atmos_phy_sf_sflx_sh(ia,ja) )
197  allocate( atmos_phy_sf_sflx_lh(ia,ja) )
198  allocate( atmos_phy_sf_sflx_gh(ia,ja) )
199  allocate( atmos_phy_sf_sflx_qtrc(ia,ja,qa) )
200  atmos_phy_sf_sflx_mw(:,:) = undef
201  atmos_phy_sf_sflx_mu(:,:) = undef
202  atmos_phy_sf_sflx_mv(:,:) = undef
203  atmos_phy_sf_sflx_sh(:,:) = undef
204  atmos_phy_sf_sflx_lh(:,:) = undef
205  atmos_phy_sf_sflx_gh(:,:) = undef
206  atmos_phy_sf_sflx_qtrc(:,:,:) = undef
207 
208  allocate( atmos_phy_sf_u10(ia,ja) )
209  allocate( atmos_phy_sf_v10(ia,ja) )
210  allocate( atmos_phy_sf_t2(ia,ja) )
211  allocate( atmos_phy_sf_q2(ia,ja) )
212  atmos_phy_sf_u10(:,:) = undef
213  atmos_phy_sf_v10(:,:) = undef
214  atmos_phy_sf_t2(:,:) = undef
215  atmos_phy_sf_q2(:,:) = undef
216 
217  !--- read namelist
218  rewind(io_fid_conf)
219  read(io_fid_conf,nml=param_atmos_phy_sf_vars,iostat=ierr)
220  if( ierr < 0 ) then !--- missing
221  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
222  elseif( ierr > 0 ) then !--- fatal error
223  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_VARS. Check!'
224  call prc_mpistop
225  endif
226  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_sf_vars)
227 
228  ! [add] 2014/08/28 A.Noda
229  atmos_phy_sf_sfc_temp(:,:) = atmos_phy_sf_default_sfc_temp
230  atmos_phy_sf_sfc_albedo(:,:,:) = atmos_phy_sf_default_sfc_albedo
231  atmos_phy_sf_sfc_z0m(:,:) = atmos_phy_sf_default_sfc_z0
232  atmos_phy_sf_sfc_z0h(:,:) = atmos_phy_sf_default_sfc_z0
233  atmos_phy_sf_sfc_z0e(:,:) = atmos_phy_sf_default_sfc_z0
234 
235  if( io_l ) write(io_fid_log,*)
236  if( io_l ) write(io_fid_log,*) '*** [ATMOS_PHY_SF] prognostic/diagnostic variables'
237  if( io_l ) write(io_fid_log,'(1x,A,A15,A,A32,3(A))') &
238  '*** |','VARNAME ','|', 'DESCRIPTION ','[', 'UNIT ',']'
239  do iv = 1, vmax
240  if( io_l ) write(io_fid_log,'(1x,A,i3,A,A15,A,A32,3(A))') &
241  '*** NO.',iv,'|',var_name(iv),'|',var_desc(iv),'[',var_unit(iv),']'
242  enddo
243 
244  if( io_l ) write(io_fid_log,*)
245  if ( atmos_phy_sf_restart_in_basename /= '' ) then
246  if( io_l ) write(io_fid_log,*) '*** Restart input? : ', trim(atmos_phy_sf_restart_in_basename)
247  else
248  if( io_l ) write(io_fid_log,*) '*** Restart input? : NO'
249  endif
251  .AND. atmos_phy_sf_restart_out_basename /= '' ) then
252  if( io_l ) write(io_fid_log,*) '*** Restart output? : ', trim(atmos_phy_sf_restart_out_basename)
253  else
254  if( io_l ) write(io_fid_log,*) '*** Restart output? : NO'
256  endif
257 
258  return
259  end subroutine atmos_phy_sf_vars_setup
260 
261  !-----------------------------------------------------------------------------
263  subroutine atmos_phy_sf_vars_fillhalo
264  use scale_const, only: &
265  i_sw => const_i_sw, &
266  i_lw => const_i_lw
267  use scale_comm, only: &
268  comm_vars8, &
269  comm_wait
270  implicit none
271  !---------------------------------------------------------------------------
272 
273  call comm_vars8( atmos_phy_sf_sfc_temp(:,:), 1 )
274  call comm_vars8( atmos_phy_sf_sfc_albedo(:,:,i_lw), 2 )
275  call comm_vars8( atmos_phy_sf_sfc_albedo(:,:,i_sw), 3 )
276  call comm_vars8( atmos_phy_sf_sfc_z0m(:,:), 4 )
277  call comm_vars8( atmos_phy_sf_sfc_z0h(:,:), 5 )
278  call comm_vars8( atmos_phy_sf_sfc_z0e(:,:), 6 )
279 
280  call comm_wait ( atmos_phy_sf_sfc_temp(:,:), 1 )
281  call comm_wait ( atmos_phy_sf_sfc_albedo(:,:,i_lw), 2 )
282  call comm_wait ( atmos_phy_sf_sfc_albedo(:,:,i_sw), 3 )
283  call comm_wait ( atmos_phy_sf_sfc_z0m(:,:), 4 )
284  call comm_wait ( atmos_phy_sf_sfc_z0h(:,:), 5 )
285  call comm_wait ( atmos_phy_sf_sfc_z0e(:,:), 6 )
286 
287  return
288  end subroutine atmos_phy_sf_vars_fillhalo
289 
290  !-----------------------------------------------------------------------------
293  use scale_const, only: &
294  i_sw => const_i_sw, &
295  i_lw => const_i_lw
296  use scale_fileio, only: &
297  fileio_read
298  use scale_rm_statistics, only: &
299  stat_total
300  implicit none
301 
302  real(RP) :: total
303  !---------------------------------------------------------------------------
304 
305  if( io_l ) write(io_fid_log,*)
306  if( io_l ) write(io_fid_log,*) '*** Input restart file (ATMOS_PHY_SF) ***'
307 
308  if ( atmos_phy_sf_restart_in_basename /= '' ) then
309  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(atmos_phy_sf_restart_in_basename)
310 
311  call fileio_read( atmos_phy_sf_sfc_temp(:,:), & ! [OUT]
312  atmos_phy_sf_restart_in_basename, var_name(1), 'XY', step=1 ) ! [IN]
313  call fileio_read( atmos_phy_sf_sfc_albedo(:,:,i_lw), & ! [OUT]
314  atmos_phy_sf_restart_in_basename, var_name(2), 'XY', step=1 ) ! [IN]
315  call fileio_read( atmos_phy_sf_sfc_albedo(:,:,i_sw), & ! [OUT]
316  atmos_phy_sf_restart_in_basename, var_name(3), 'XY', step=1 ) ! [IN]
317  call fileio_read( atmos_phy_sf_sfc_z0m(:,:), & ! [OUT]
318  atmos_phy_sf_restart_in_basename, var_name(4), 'XY', step=1 ) ! [IN]
319  call fileio_read( atmos_phy_sf_sfc_z0h(:,:), & ! [OUT]
320  atmos_phy_sf_restart_in_basename, var_name(5), 'XY', step=1 ) ! [IN]
321  call fileio_read( atmos_phy_sf_sfc_z0e(:,:), & ! [OUT]
322  atmos_phy_sf_restart_in_basename, var_name(6), 'XY', step=1 ) ! [IN]
323 
325 
326  call stat_total( total, atmos_phy_sf_sfc_temp(:,:), var_name(1) )
327  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_lw), var_name(2) )
328  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_sw), var_name(3) )
329  call stat_total( total, atmos_phy_sf_sfc_z0m(:,:), var_name(4) )
330  call stat_total( total, atmos_phy_sf_sfc_z0h(:,:), var_name(5) )
331  call stat_total( total, atmos_phy_sf_sfc_z0e(:,:), var_name(6) )
332  else
333  if( io_l ) write(io_fid_log,*) '*** restart file for ATMOS_PHY_SF is not specified.'
334  atmos_phy_sf_sfc_temp(:,:) = atmos_phy_sf_default_sfc_temp
335  atmos_phy_sf_sfc_albedo(:,:,:) = atmos_phy_sf_default_sfc_albedo
336  atmos_phy_sf_sfc_z0m(:,:) = atmos_phy_sf_default_sfc_z0
337  atmos_phy_sf_sfc_z0h(:,:) = atmos_phy_sf_default_sfc_z0
338  atmos_phy_sf_sfc_z0e(:,:) = atmos_phy_sf_default_sfc_z0
339  endif
340 
341  return
342  end subroutine atmos_phy_sf_vars_restart_read
343 
344  !-----------------------------------------------------------------------------
347  use scale_const, only: &
348  i_sw => const_i_sw, &
349  i_lw => const_i_lw
350  use scale_time, only: &
352  use scale_fileio, only: &
353  fileio_write
354  use scale_rm_statistics, only: &
355  stat_total
356  implicit none
357 
358  character(len=20) :: timelabel
359  character(len=H_LONG) :: basename
360 
361  real(RP) :: total
362  !---------------------------------------------------------------------------
363 
364  if ( atmos_phy_sf_restart_out_basename /= '' ) then
365 
366  call time_gettimelabel( timelabel )
367  write(basename,'(A,A,A)') trim(atmos_phy_sf_restart_out_basename), '_', trim(timelabel)
368 
369  if( io_l ) write(io_fid_log,*)
370  if( io_l ) write(io_fid_log,*) '*** Output restart file (ATMOS_PHY_SF) ***'
371  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(basename)
372 
374 
375  call stat_total( total, atmos_phy_sf_sfc_temp(:,:), var_name(1) )
376  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_lw), var_name(2) )
377  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_sw), var_name(3) )
378  call stat_total( total, atmos_phy_sf_sfc_z0m(:,:), var_name(4) )
379  call stat_total( total, atmos_phy_sf_sfc_z0h(:,:), var_name(5) )
380  call stat_total( total, atmos_phy_sf_sfc_z0e(:,:), var_name(6) )
381 
382  call fileio_write( atmos_phy_sf_sfc_temp(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
383  var_name(1), var_desc(1), var_unit(1), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
384  call fileio_write( atmos_phy_sf_sfc_albedo(:,:,i_lw), basename, atmos_phy_sf_restart_out_title, & ! [IN]
385  var_name(2), var_desc(2), var_unit(2), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
386  call fileio_write( atmos_phy_sf_sfc_albedo(:,:,i_sw), basename, atmos_phy_sf_restart_out_title, & ! [IN]
387  var_name(3), var_desc(3), var_unit(3), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
388  call fileio_write( atmos_phy_sf_sfc_z0m(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
389  var_name(4), var_desc(4), var_unit(4), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
390  call fileio_write( atmos_phy_sf_sfc_z0h(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
391  var_name(5), var_desc(5), var_unit(5), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
392  call fileio_write( atmos_phy_sf_sfc_z0e(:,:), basename, atmos_phy_sf_restart_out_title, & ! [IN]
393  var_name(6), var_desc(6), var_unit(6), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
394  endif
395 
396  return
397  end subroutine atmos_phy_sf_vars_restart_write
398 
399  !-----------------------------------------------------------------------------
401  subroutine atmos_phy_sf_vars_external_in( &
402  sfc_temp_in, &
403  albedo_lw_in, &
404  albedo_sw_in, &
405  sfc_z0m_in, &
406  sfc_z0h_in, &
407  sfc_z0e_in )
408  use scale_const, only: &
409  i_sw => const_i_sw, &
410  i_lw => const_i_lw
411  implicit none
412 
413  real(RP), intent(in) :: sfc_temp_in (:,:)
414  real(RP), intent(in) :: albedo_lw_in(:,:)
415  real(RP), intent(in) :: albedo_sw_in(:,:)
416  real(RP), intent(in) :: sfc_z0m_in (:,:)
417  real(RP), intent(in) :: sfc_z0h_in (:,:)
418  real(RP), intent(in) :: sfc_z0e_in (:,:)
419  !---------------------------------------------------------------------------
420 
421  if( io_l ) write(io_fid_log,*)
422  if( io_l ) write(io_fid_log,*) '*** External Input (PHY_SF) ***'
423 
424  atmos_phy_sf_sfc_temp(:,:) = sfc_temp_in(:,:)
425  atmos_phy_sf_sfc_albedo(:,:,i_lw) = albedo_lw_in(:,:)
426  atmos_phy_sf_sfc_albedo(:,:,i_sw) = albedo_sw_in(:,:)
427  atmos_phy_sf_sfc_z0m(:,:) = sfc_z0m_in(:,:)
428  atmos_phy_sf_sfc_z0h(:,:) = sfc_z0h_in(:,:)
429  atmos_phy_sf_sfc_z0e(:,:) = sfc_z0e_in(:,:)
430 
432 
433  return
434  end subroutine atmos_phy_sf_vars_external_in
435 
436  !-----------------------------------------------------------------------------
439  use scale_const, only: &
440  i_sw => const_i_sw, &
441  i_lw => const_i_lw
442  use scale_time, only: &
444  use scale_fileio, only: &
446  implicit none
447 
448  character(len=20) :: timelabel
449  character(len=H_LONG) :: basename
450 
451  !---------------------------------------------------------------------------
452 
453  if ( atmos_phy_sf_restart_out_basename /= '' ) then
454 
455  call time_gettimelabel( timelabel )
456  write(basename,'(A,A,A)') trim(atmos_phy_sf_restart_out_basename), '_', trim(timelabel)
457 
458  if( io_l ) write(io_fid_log,*)
459  if( io_l ) write(io_fid_log,*) '*** Output restart file (ATMOS_PHY_SF) ***'
460  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(basename)
461 
462  call fileio_create(restart_fid, basename, atmos_phy_sf_restart_out_title, & ! [IN]
464  endif
465 
466  return
467  end subroutine atmos_phy_sf_vars_restart_create
468 
469  !-----------------------------------------------------------------------------
472  use scale_fileio, only: &
474  implicit none
475 
476  if ( restart_fid .NE. -1 ) then
477  call fileio_enddef( restart_fid ) ! [IN]
478  endif
479 
480  return
481  end subroutine atmos_phy_sf_vars_restart_enddef
482 
483  !-----------------------------------------------------------------------------
486  use scale_fileio, only: &
488  implicit none
489 
490  if ( restart_fid .NE. -1 ) then
491  call fileio_close( restart_fid ) ! [IN]
492  restart_fid = -1
493  endif
494 
495  return
496  end subroutine atmos_phy_sf_vars_restart_close
497 
498  !-----------------------------------------------------------------------------
501  use scale_fileio, only: &
503  implicit none
504 
505  !---------------------------------------------------------------------------
506 
507  if ( restart_fid .NE. -1 ) then
508 
509  call fileio_def_var( restart_fid, var_id(1), var_name(1), var_desc(1), &
510  var_unit(1), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
511  call fileio_def_var( restart_fid, var_id(2), var_name(2), var_desc(2), &
512  var_unit(2), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
513  call fileio_def_var( restart_fid, var_id(3), var_name(3), var_desc(3), &
514  var_unit(3), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
515  call fileio_def_var( restart_fid, var_id(4), var_name(4), var_desc(4), &
516  var_unit(4), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
517  call fileio_def_var( restart_fid, var_id(5), var_name(5), var_desc(5), &
518  var_unit(5), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
519  call fileio_def_var( restart_fid, var_id(6), var_name(6), var_desc(6), &
520  var_unit(6), 'XY', atmos_phy_sf_restart_out_dtype ) ! [IN]
521  endif
522 
523  return
524  end subroutine atmos_phy_sf_vars_restart_def_var
525 
526  !-----------------------------------------------------------------------------
529  use scale_const, only: &
530  i_sw => const_i_sw, &
531  i_lw => const_i_lw
532  use scale_fileio, only: &
533  fileio_write_var
534  use scale_rm_statistics, only: &
535  stat_total
536  implicit none
537 
538  real(RP) :: total
539  !---------------------------------------------------------------------------
540 
541  if ( restart_fid .NE. -1 ) then
542 
544 
545  call stat_total( total, atmos_phy_sf_sfc_temp(:,:), var_name(1) )
546  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_lw), var_name(2) )
547  call stat_total( total, atmos_phy_sf_sfc_albedo(:,:,i_sw), var_name(3) )
548  call stat_total( total, atmos_phy_sf_sfc_z0m(:,:), var_name(4) )
549  call stat_total( total, atmos_phy_sf_sfc_z0h(:,:), var_name(5) )
550  call stat_total( total, atmos_phy_sf_sfc_z0e(:,:), var_name(6) )
551 
552  call fileio_write_var( restart_fid, var_id(1), atmos_phy_sf_sfc_temp(:,:), & ! [IN]
553  var_name(1), 'XY' ) ! [IN]
554  call fileio_write_var( restart_fid, var_id(2), atmos_phy_sf_sfc_albedo(:,:,i_lw), & ! [IN]
555  var_name(2), 'XY' ) ! [IN]
556  call fileio_write_var( restart_fid, var_id(3), atmos_phy_sf_sfc_albedo(:,:,i_sw), & ! [IN]
557  var_name(3), 'XY' ) ! [IN]
558  call fileio_write_var( restart_fid, var_id(4), atmos_phy_sf_sfc_z0m(:,:), & ! [IN]
559  var_name(4), 'XY' ) ! [IN]
560  call fileio_write_var( restart_fid, var_id(5), atmos_phy_sf_sfc_z0h(:,:), & ! [IN]
561  var_name(5), 'XY' ) ! [IN]
562  call fileio_write_var( restart_fid, var_id(6), atmos_phy_sf_sfc_z0e(:,:), & ! [IN]
563  var_name(6), 'XY' ) ! [IN]
564  endif
565 
566  return
568 
569 end module mod_atmos_phy_sf_vars
subroutine, public atmos_phy_sf_vars_restart_def_var
Write restart.
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
subroutine, public atmos_phy_sf_vars_restart_enddef
Exit netCDF define mode.
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
subroutine, public atmos_phy_sf_vars_restart_write_var
Write variables to restart file.
module grid index
module ATMOSPHERIC Surface Variables
subroutine, public atmos_phy_sf_vars_restart_write
Write restart.
subroutine, public atmos_phy_sf_vars_restart_close
Close restart file.
module TRACER
subroutine, public fileio_def_var(fid, vid, varname, desc, unit, axistype, datatype, timeintv)
Define a variable to file.
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
subroutine, public fileio_create(fid, basename, title, datatype, date, subsec, append, nozcoord)
Create/open a netCDF file.
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
subroutine, public fileio_enddef(fid)
Exit netCDF file define mode.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momy_t
subroutine, public atmos_phy_sf_vars_restart_create
Create restart file.
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
subroutine, public fileio_close(fid)
Close a netCDF file.
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?