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  logical, public :: atmos_phy_sf_restart_in_postfix_timelabel = .false.
52  character(len=H_LONG), public :: atmos_phy_sf_restart_out_basename = ''
53  logical, public :: atmos_phy_sf_restart_out_postfix_timelabel = .true.
54  character(len=H_MID), public :: atmos_phy_sf_restart_out_title = 'ATMOS_PHY_SF restart'
55  character(len=H_SHORT), public :: atmos_phy_sf_restart_out_dtype = 'DEFAULT'
56 
57  real(RP), public, allocatable :: atmos_phy_sf_dens_t(:,:) ! tendency DENS [ kg/m3/s]
58  real(RP), public, allocatable :: atmos_phy_sf_momz_t(:,:) ! tendency MOMZ [m/s*kg/m3/s]
59  real(RP), public, allocatable :: atmos_phy_sf_momx_t(:,:) ! tendency MOMX [m/s*kg/m3/s]
60  real(RP), public, allocatable :: atmos_phy_sf_momy_t(:,:) ! tendency MOMY [m/s*kg/m3/s]
61  real(RP), public, allocatable :: atmos_phy_sf_rhot_t(:,:) ! tendency RHOT [K *kg/m3/s]
62  real(RP), public, allocatable :: atmos_phy_sf_rhoq_t(:,:,:) ! tendency rho*QTRC [ kg/kg/s]
63 
64  real(RP), public, allocatable :: atmos_phy_sf_sfc_temp (:,:) ! surface skin temperature [K]
65  real(RP), public, allocatable :: atmos_phy_sf_sfc_albedo(:,:,:) ! surface albedo (0-1)
66  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0m (:,:) ! surface roughness length, ocean only [m]
67  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0h (:,:) ! surface roughness length, ocean only [m]
68  real(RP), public, allocatable :: atmos_phy_sf_sfc_z0e (:,:) ! surface roughness length, ocean only [m]
69 
70  ! surface diagnostic variables
71  real(RP), public, allocatable :: atmos_phy_sf_sfc_dens (:,:) ! surface atmosphere density [kg/m3]
72  real(RP), public, allocatable :: atmos_phy_sf_sfc_pres (:,:) ! surface atmosphere pressure [Pa]
73 
74  real(RP), public, allocatable :: atmos_phy_sf_sflx_mw (:,:) ! z-momentum flux (area center) [m/s*kg/m2/s]
75  real(RP), public, allocatable :: atmos_phy_sf_sflx_mu (:,:) ! x-momentum flux (area center) [m/s*kg/m2/s]
76  real(RP), public, allocatable :: atmos_phy_sf_sflx_mv (:,:) ! y-momentum flux (area center) [m/s*kg/m2/s]
77  real(RP), public, allocatable :: atmos_phy_sf_sflx_sh (:,:) ! sensible heat flux [J/m2/s]
78  real(RP), public, allocatable :: atmos_phy_sf_sflx_lh (:,:) ! latent heat flux [J/m2/s]
79  real(RP), public, allocatable :: atmos_phy_sf_sflx_gh (:,:) ! ground heat flux [J/m2/s]
80  real(RP), public, allocatable :: atmos_phy_sf_sflx_qtrc (:,:,:) ! tracer mass flux [kg/m2/s]
81 
82  real(RP), public, allocatable :: atmos_phy_sf_u10 (:,:) ! 10m x-wind [m/s]
83  real(RP), public, allocatable :: atmos_phy_sf_v10 (:,:) ! 10m y-wind [m/s]
84  real(RP), public, allocatable :: atmos_phy_sf_t2 (:,:) ! 2m temperature [K]
85  real(RP), public, allocatable :: atmos_phy_sf_q2 (:,:) ! 2m specific humidity [kg/kg]
86 
87 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_QEMIS(:,:,:) ! tracer emission flux [kg/m2/s]
88 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_QDEP (:,:,:) ! tracer deposition flux [kg/m2/s]
89 ! real(RP), public, allocatable :: ATMOS_PHY_SF_SFLX_VDEP (:,:,:) ! tracer deposition velocity [m/s]
90 
91  !-----------------------------------------------------------------------------
92  !
93  !++ Private procedure
94  !
95  !-----------------------------------------------------------------------------
96  !
97  !++ Private parameters & variables
98  !
99  integer, private, parameter :: vmax = 6
100  integer, private, parameter :: i_sfc_temp = 1
101  integer, private, parameter :: i_sfc_albedo_lw = 2
102  integer, private, parameter :: i_sfc_albedo_sw = 3
103  integer, private, parameter :: i_sfc_z0m = 4
104  integer, private, parameter :: i_sfc_z0h = 5
105  integer, private, parameter :: i_sfc_z0e = 6
106 
107  character(len=H_SHORT), private :: var_name(vmax)
108  character(len=H_MID), private :: var_desc(vmax)
109  character(len=H_SHORT), private :: var_unit(vmax)
110  integer, private :: var_id(vmax)
111  integer, private :: restart_fid = -1 ! file ID
112 
113  data var_name / 'SFC_TEMP', &
114  'SFC_ALB_LW', &
115  'SFC_ALB_SW', &
116  'SFC_Z0M', &
117  'SFC_Z0H', &
118  'SFC_Z0E' /
119 
120  data var_desc / 'surface skin temperature', &
121  'surface albedo for longwave', &
122  'surface albedo for shortwave', &
123  'surface roughness length (momentum)', &
124  'surface roughness length (heat)', &
125  'surface roughness length (vapor)' /
126 
127  data var_unit / 'K', &
128  '1', &
129  '1', &
130  'm', &
131  'm', &
132  'm' /
133 
134  real(RP), private :: atmos_phy_sf_default_sfc_temp = 300.0_rp
135  real(RP), private :: atmos_phy_sf_default_sfc_albedo = 0.4_rp
136  real(RP), private :: atmos_phy_sf_default_sfc_z0 = 1e-4_rp
137 
138  !-----------------------------------------------------------------------------
139 contains
140  !-----------------------------------------------------------------------------
142  subroutine atmos_phy_sf_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 / &
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
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'
263  endif
264 
265  return
266  end subroutine atmos_phy_sf_vars_setup
267 
268  !-----------------------------------------------------------------------------
270  subroutine atmos_phy_sf_vars_fillhalo
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
295  end subroutine atmos_phy_sf_vars_fillhalo
296 
297  !-----------------------------------------------------------------------------
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 
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
336  end subroutine atmos_phy_sf_vars_restart_open
337 
338  !-----------------------------------------------------------------------------
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
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
391  end subroutine atmos_phy_sf_vars_restart_read
392 
393  !-----------------------------------------------------------------------------
395  subroutine atmos_phy_sf_vars_external_in( &
396  sfc_temp_in, &
397  albedo_lw_in, &
398  albedo_sw_in, &
399  sfc_z0m_in, &
400  sfc_z0h_in, &
401  sfc_z0e_in )
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 
426 
427  return
428  end subroutine atmos_phy_sf_vars_external_in
429 
430  !-----------------------------------------------------------------------------
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 
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]
459  endif
460 
461  return
462  end subroutine atmos_phy_sf_vars_restart_create
463 
464  !-----------------------------------------------------------------------------
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
476  end subroutine atmos_phy_sf_vars_restart_enddef
477 
478  !-----------------------------------------------------------------------------
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
496  end subroutine atmos_phy_sf_vars_restart_close
497 
498  !-----------------------------------------------------------------------------
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
523  end subroutine atmos_phy_sf_vars_restart_def_var
524 
525  !-----------------------------------------------------------------------------
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 
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
570  end subroutine atmos_phy_sf_vars_restart_write
571 
572 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_open
Open restart file for read.
logical, public statistics_checktotal
calc&report variable totals to logfile?
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:95
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:61
subroutine, public fileio_flush(fid)
Flush all pending requests to a netCDF file (PnetCDF only)
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.
real(rp), public const_undef
Definition: scale_const.F90:43
module Statistics
module grid index
module ATMOSPHERIC Surface Variables
subroutine, public atmos_phy_sf_vars_restart_write
Write variables to restart file.
subroutine, public atmos_phy_sf_vars_restart_close
Close restart file.
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:62
module TRACER
integer, public ia
of whole cells: x, 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 input file.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
module COMMUNICATION
Definition: scale_comm.F90:23
character(len=h_short), public atmos_phy_sf_restart_out_dtype
REAL4 or REAL8.
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
logical, public atmos_phy_sf_restart_out_postfix_timelabel
Add timelabel to the basename of output file?
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_aggregate
do parallel I/O through PnetCDF
Definition: scale_stdio.F90:66
subroutine, public fileio_open(fid, basename)
open a netCDF file for read
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:96
module PRECISION
logical, public atmos_phy_sf_restart_in_postfix_timelabel
Add timelabel to the basename of input file?
subroutine, public fileio_def_var(fid, vid, varname, desc, unit, axistype, datatype, timeintv, nsteps)
Define a variable to file.
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 io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
integer, public ja
of whole cells: y, local, with HALO
logical, public atmos_phy_sf_restart_output
output restart file?