SCALE-RM
scale_history.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
15 !-------------------------------------------------------------------------------
17  !-----------------------------------------------------------------------------
18  !
19  !++ Used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
27  !-----------------------------------------------------------------------------
28  implicit none
29  private
30  !-----------------------------------------------------------------------------
31  !
32  !++ Public procedures
33  !
34  public :: hist_setup
35  public :: hist_reg
36  public :: hist_query
37  public :: hist_in
38  public :: hist_put
39  public :: hist_get
40  public :: hist_write
41  public :: hist_switch
42 
43  interface hist_in
44  module procedure hist_in_0d
45  module procedure hist_in_1d
46  module procedure hist_in_2d
47  module procedure hist_in_3d
48  end interface hist_in
49  interface hist_put
50  module procedure hist_put_0d
51  module procedure hist_put_1d
52  module procedure hist_put_2d
53  module procedure hist_put_3d
54  end interface hist_put
55  interface hist_get
56  module procedure hist_get_1d
57  module procedure hist_get_2d
58  module procedure hist_get_3d
59  end interface hist_get
60 
61  !-----------------------------------------------------------------------------
62  !
63  !++ Public parameters & variables
64  !
65  !-----------------------------------------------------------------------------
66  !
67  !++ Private procedures
68  !
69  private :: hist_put_axes
70 
71  !-----------------------------------------------------------------------------
72  !
73  !++ Private parameters & variables
74  !
75  logical :: enabled
76  integer :: im, jm
77  integer :: ims, ime
78  integer :: jms, jme
79  !-----------------------------------------------------------------------------
80 contains
81  !-----------------------------------------------------------------------------
83  subroutine hist_setup
84  use gtool_history, only: &
86  use scale_process, only: &
87  prc_mpistop, &
90  use scale_rm_process, only: &
92  use scale_time, only: &
93  time_nowdate, &
94  time_nowms, &
95  time_dtsec, &
98  implicit none
99 
100  character(len=H_MID) :: HISTORY_H_TITLE = 'SCALE-RM HISTORY OUTPUT'
101  character(len=H_SHORT) :: HISTORY_T_UNITS = 'seconds'
102  character(len=H_MID) :: HISTORY_T_SINCE = ''
103 
104  logical :: HIST_BND = .true.
105 
106  namelist / param_hist / &
107  hist_bnd
108 
109  real(DP) :: start_daysec
110  integer :: rankidx(2)
111  integer :: ierr
112  !---------------------------------------------------------------------------
113 
114  if( io_l ) write(io_fid_log,*)
115  if( io_l ) write(io_fid_log,*) '++++++ Module[HISTORY] / Categ[ATMOS-RM IO] / Origin[SCALElib]'
116 
117  !--- read namelist
118  rewind(io_fid_conf)
119  read(io_fid_conf,nml=param_hist,iostat=ierr)
120  if( ierr < 0 ) then !--- missing
121  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
122  elseif( ierr > 0 ) then !--- fatal error
123  write(*,*) 'xxx Not appropriate names in namelist PARAM_HIST. Check!'
124  call prc_mpistop
125  endif
126  if( io_lnml ) write(io_fid_log,nml=param_hist)
127 
128 
129  call prof_rapstart('FILE_O_NetCDF', 2)
130 
131  rankidx(1) = prc_2drank(prc_myrank, 1)
132  rankidx(2) = prc_2drank(prc_myrank, 2)
133 
134  start_daysec = time_startdaysec
135 ! if ( TIME_OFFSET_YEAR > 0 ) then
136 ! HISTORY_T_SINCE = '0000-01-01 00:00:00'
137 ! write(HISTORY_T_SINCE(1:4),'(i4)') TIME_OFFSET_YEAR
138 ! end if
139  if ( time_nowdate(1) > 0 ) then
140  write(history_t_since, '(i4.4,a1,i2.2,a1,i2.2,a1,i2.2,a1,i2.2,a1,i2.2)') &
141  time_nowdate(1),'-',time_nowdate(2),'-',time_nowdate(3),' ', &
142  time_nowdate(4),':',time_nowdate(5),':',time_nowdate(6)
143  start_daysec = time_nowms
144  end if
145 
146  if ( hist_bnd ) then
147  im = imaxb
148  jm = jmaxb
149  ims = isb
150  ime = ieb
151  jms = jsb
152  jme = jeb
153  else
154  im = imax
155  jm = jmax
156  ims = is
157  ime = ie
158  jms = js
159  jme = je
160  end if
161 
162  call historyinit( history_h_title, &
163  h_source, &
164  h_institute, &
165  im*jm*kmax, &
166  prc_masterrank, &
167  prc_myrank, &
168  rankidx, &
169  start_daysec, &
170  time_dtsec, &
171  time_units = history_t_units, &
172  time_since = history_t_since, &
173  namelist_fid = io_fid_conf )
174 
175  call hist_put_axes
176 
177  call prof_rapend ('FILE_O_NetCDF', 2)
178 
179  enabled = .true.
180 
181  return
182  end subroutine hist_setup
183 
184  !-----------------------------------------------------------------------------
186  subroutine hist_switch( switch )
187  logical, intent(in) :: switch
188 
189  enabled = switch
190 
191  return
192  end subroutine hist_switch
193 
194  !-----------------------------------------------------------------------------
196  subroutine hist_put_axes
197  use gtool_history, only: &
198  historyputaxis, &
199  historyputassociatedcoordinates
200  use scale_const, only: &
201  d2r => const_d2r
202  use scale_grid, only: &
203  grid_cz, &
204  grid_cx, &
205  grid_cy, &
206  grid_fz, &
207  grid_fx, &
208  grid_fy, &
209  grid_cdz, &
210  grid_cdx, &
211  grid_cdy, &
212  grid_fdz, &
213  grid_fdx, &
214  grid_fdy, &
215  grid_cbfz, &
216  grid_cbfx, &
217  grid_cbfy, &
218  grid_fbfz, &
219  grid_fbfx, &
220  grid_fbfy, &
221  grid_cxg, &
222  grid_cyg, &
223  grid_fxg, &
224  grid_fyg, &
225  grid_cbfxg, &
226  grid_cbfyg, &
227  grid_fbfxg, &
228  grid_fbfyg
229  use scale_land_grid, only: &
230  grid_lcz, &
231  grid_lfz, &
232  grid_lcdz
233  use scale_urban_grid, only: &
234  grid_ucz, &
235  grid_ufz, &
236  grid_ucdz
237  use scale_grid_real, only: &
238  real_cz, &
239  real_fz, &
240  real_lon, &
241  real_lonx, &
242  real_lony, &
243  real_lonxy, &
244  real_lat, &
245  real_latx, &
246  real_laty, &
247  real_latxy
248  use scale_topography, only: &
249  topo_zsfc
250  use scale_landuse, only: &
252  implicit none
253 
254  real(RP) :: AXIS (im,jm,kmax)
255  character(len=2) :: AXIS_name(3)
256 
257  integer :: k, i, j
258  !---------------------------------------------------------------------------
259 
260  call historyputaxis( 'x', 'X', 'm', 'x', grid_cx(ims:ime) )
261  call historyputaxis( 'y', 'Y', 'm', 'y', grid_cy(jms:jme) )
262  call historyputaxis( 'z', 'Z', 'm', 'z', grid_cz(ks:ke) )
263  call historyputaxis( 'xh', 'X (half level)', 'm', 'xh', grid_fx(ims:ime) )
264  call historyputaxis( 'yh', 'Y (half level)', 'm', 'yh', grid_fy(jms:jme) )
265  call historyputaxis( 'zh', 'Z (half level)', 'm', 'zh', grid_fz(ks:ke) )
266 
267  call historyputaxis( 'lz', 'LZ', 'm', 'lz', grid_lcz(lks:lke), down=.true. )
268  call historyputaxis( 'lzh', 'LZ (half level)', 'm', 'lzh', grid_lfz(lks:lke), down=.true. )
269 
270  call historyputaxis( 'uz', 'UZ', 'm', 'uz', grid_ucz(uks:uke), down=.true. )
271  call historyputaxis( 'uzh', 'UZ (half level)', 'm', 'uzh', grid_ufz(uks:uke), down=.true. )
272 
273  call historyputaxis( 'CZ', 'Atmos Grid Center Position Z', 'm', 'CZ', grid_cz )
274  call historyputaxis( 'CX', 'Atmos Grid Center Position X', 'm', 'CX', grid_cx )
275  call historyputaxis( 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', grid_cy )
276  call historyputaxis( 'FZ', 'Atmos Grid Face Position Z', 'm', 'FZ', grid_fz )
277  call historyputaxis( 'FX', 'Atmos Grid Face Position X', 'm', 'FX', grid_fx )
278  call historyputaxis( 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', grid_fy )
279 
280  call historyputaxis( 'CDZ', 'Grid Cell length Z', 'm', 'CZ', grid_cdz )
281  call historyputaxis( 'CDX', 'Grid Cell length X', 'm', 'CX', grid_cdx )
282  call historyputaxis( 'CDY', 'Grid Cell length Y', 'm', 'CY', grid_cdy )
283  call historyputaxis( 'FDZ', 'Grid distance Z', 'm', 'FDZ', grid_fdz )
284  call historyputaxis( 'FDX', 'Grid distance X', 'm', 'FDX', grid_fdx )
285  call historyputaxis( 'FDY', 'Grid distance Y', 'm', 'FDY', grid_fdy )
286 
287  call historyputaxis( 'LCZ', 'Land Grid Center Position Z', 'm', 'LCZ', grid_lcz, down=.true. )
288  call historyputaxis( 'LFZ', 'Land Grid Face Position Z', 'm', 'LFZ', grid_lfz, down=.true. )
289  call historyputaxis( 'LCDZ', 'Land Grid Cell length Z', 'm', 'LCZ', grid_lcdz )
290 
291  call historyputaxis( 'UCZ', 'Urban Grid Center Position Z', 'm', 'UCZ', grid_ucz, down=.true. )
292  call historyputaxis( 'UFZ', 'Urban Grid Face Position Z', 'm', 'UFZ', grid_ufz , down=.true. )
293  call historyputaxis( 'UCDZ', 'Urban Grid Cell length Z', 'm', 'UCZ', grid_ucdz )
294 
295  call historyputaxis('CBFZ', 'Boundary factor Center Z', '1', 'CZ', grid_cbfz)
296  call historyputaxis('CBFX', 'Boundary factor Center X', '1', 'CX', grid_cbfx)
297  call historyputaxis('CBFY', 'Boundary factor Center Y', '1', 'CY', grid_cbfy)
298  call historyputaxis('FBFZ', 'Boundary factor Face Z', '1', 'CZ', grid_fbfz)
299  call historyputaxis('FBFX', 'Boundary factor Face X', '1', 'CX', grid_fbfx)
300  call historyputaxis('FBFY', 'Boundary factor Face Y', '1', 'CY', grid_fbfy)
301 
302  call historyputaxis('CXG', 'Grid Center Position X (global)', 'm', 'CXG', grid_cxg)
303  call historyputaxis('CYG', 'Grid Center Position Y (global)', 'm', 'CYG', grid_cyg)
304  call historyputaxis('FXG', 'Grid Face Position X (global)', 'm', 'FXG', grid_fxg)
305  call historyputaxis('FYG', 'Grid Face Position Y (global)', 'm', 'FYG', grid_fyg)
306 
307  call historyputaxis('CBFXG', 'Boundary factor Center X (global)', '1', 'CXG', grid_cbfxg)
308  call historyputaxis('CBFYG', 'Boundary factor Center Y (global)', '1', 'CYG', grid_cbfyg)
309  call historyputaxis('FBFXG', 'Boundary factor Face X (global)', '1', 'CXG', grid_fbfxg)
310  call historyputaxis('FBFYG', 'Boundary factor Face Y (global)', '1', 'CYG', grid_fbfyg)
311 
312  ! associate coordinates
313  do k = 1, kmax
314  axis(1:im,1:jm,k) = real_cz(k+ks-1,ims:ime,jms:jme)
315  enddo
316  axis_name(1:3) = (/'x ','y ', 'z '/)
317  call historyputassociatedcoordinates( 'height' , 'height' , &
318  'm' , axis_name(1:3), axis(:,:,:) )
319 
320  do k = 1, kmax
321  axis(1:im,1:jm,k) = real_fz(k+ks-1,ims:ime,jms:jme)
322  enddo
323  axis_name(1:3) = (/'x ','y ', 'zh'/)
324  call historyputassociatedcoordinates( 'height_xyw' , 'height (half level xyw)' , &
325  'm' , axis_name(1:3), axis(:,:,:) )
326 
327  do k = 1, kmax
328  do j = 1, jm
329  do i = 1, im-1
330  axis(i,j,k) = ( real_cz(k+ks-1,ims+i-1,jms+j-1) + real_cz(k+ks-1,ims+i,jms+j-1) ) * 0.5_rp
331  enddo
332  enddo
333  enddo
334  do k = 1, kmax
335  do j = 1, jm
336  axis(im,j,k) = real_cz(k+ks-1,ims+im-1,jms+j-1)
337  enddo
338  enddo
339  axis_name(1:3) = (/'xh','y ', 'z '/)
340  call historyputassociatedcoordinates( 'height_uyz' , 'height (half level uyz)' , &
341  'm' , axis_name(1:3), axis(:,:,:) )
342 
343  do k = 1, kmax
344  do j = 1, jm-1
345  do i = 1, im
346  axis(i,j,k) = ( real_cz(k+ks-1,ims+i-1,jms+j-1) + real_cz(k+ks-1,ims+i-1,jms+j) ) * 0.5_rp
347  enddo
348  enddo
349  enddo
350  do k = 1, kmax
351  do i = 1, im
352  axis(i,jm,k) = real_cz(k+ks-1,ims+i-1,jms+jm-1)
353  enddo
354  enddo
355  axis_name(1:3) = (/'x ','yh', 'z '/)
356  call historyputassociatedcoordinates( 'height_xvz' , 'height (half level xvz)' , &
357  'm' , axis_name(1:3), axis(:,:,:) )
358 
359  do k = 1, kmax
360  do j = 1, jm-1
361  do i = 1, im-1
362  axis(i,j,k) = ( real_cz(k+ks-1,ims+i-1,jms+j-1) + real_cz(k+ks-1,ims+i ,jms+j-1) &
363  + real_cz(k+ks-1,ims+i-1,jms+j ) + real_cz(k+ks-1,ims+i ,jms+j ) ) * 0.25_rp
364  enddo
365  enddo
366  enddo
367  do k = 1, kmax
368  do i = 1, im
369  axis(i,jm,k) = real_cz(k+ks-1,ims+i-1,jms+jm-1)
370  enddo
371  enddo
372  do k = 1, kmax
373  do j = 1, jm
374  axis(im,j,k) = real_cz(k+ks-1,ims+im-1,jms+j-1)
375  enddo
376  enddo
377  axis_name(1:3) = (/'xh','yh', 'z '/)
378  call historyputassociatedcoordinates( 'height_uvz' , 'height (half level uvz)' , &
379  'm' , axis_name(1:3), axis(:,:,:) )
380 
381  do k = 1, kmax
382  do j = 1, jm
383  do i = 1, im-1
384  axis(i,j,k) = ( real_fz(k+ks-1,ims+i-1,jms+j-1) + real_fz(k+ks-1,ims+i,jms+j-1) ) * 0.5_rp
385  enddo
386  enddo
387  enddo
388  do k = 1, kmax
389  do j = 1, jm
390  axis(im,j,k) = real_fz(k+ks-1,ims+im-1,jms+j-1)
391  enddo
392  enddo
393  axis_name(1:3) = (/'xh','y ', 'zh'/)
394  call historyputassociatedcoordinates( 'height_uyw' , 'height (half level uyw)' , &
395  'm' , axis_name(1:3), axis(:,:,:) )
396 
397  do k = 1, kmax
398  do j = 1, jm-1
399  do i = 1, im
400  axis(i,j,k) = ( real_fz(k+ks-1,ims+i-1,jms+j-1) + real_fz(k+ks-1,ims+i-1,jms+j) ) * 0.5_rp
401  enddo
402  enddo
403  enddo
404  do k = 1, kmax
405  do i = 1, im
406  axis(i,jm,k) = real_fz(k+ks-1,ims+i-1,jms+jm-1)
407  enddo
408  enddo
409  axis_name(1:3) = (/'x ','yh', 'zh'/)
410  call historyputassociatedcoordinates( 'height_xvw' , 'height (half level xvw)' , &
411  'm' , axis_name(1:3), axis(:,:,:) )
412 
413  do k = 1, kmax
414  do j = 1, jm-1
415  do i = 1, im-1
416  axis(i,j,k) = ( real_fz(k+ks-1,ims+i-1,jms+j-1) + real_fz(k+ks-1,ims+i ,jms+j-1) &
417  + real_fz(k+ks-1,ims+i-1,jms+j ) + real_fz(k+ks-1,ims+i ,jms+j ) ) * 0.25_rp
418  enddo
419  enddo
420  enddo
421  do k = 1, kmax
422  do i = 1, im
423  axis(i,jm,k) = real_fz(k+ks-1,ims+i-1,jms+jm-1)
424  enddo
425  enddo
426  do k = 1, kmax
427  do j = 1, jm
428  axis(im,j,k) = real_fz(k+ks-1,ims+im-1,jms+j-1)
429  enddo
430  enddo
431  axis_name(1:3) = (/'xh','yh', 'zh'/)
432  call historyputassociatedcoordinates( 'height_uvw' , 'height (half level uvw)' , &
433  'm' , axis_name(1:3), axis(:,:,:) )
434 
435  axis(1:im,1:jm,1) = real_lon(ims:ime,jms:jme) / d2r
436  axis_name(1:2) = (/'x ','y '/)
437  call historyputassociatedcoordinates( 'lon' , 'longitude' , &
438  'degrees_east' , axis_name(1:2), axis(:,:,1) )
439 
440  axis(1:im,1:jm,1) = real_lonx(ims:ime,jms:jme) / d2r
441  axis_name(1:2) = (/'xh','y '/)
442  call historyputassociatedcoordinates( 'lon_uy', 'longitude (half level uy)', &
443  'degrees_east' , axis_name(1:2), axis(:,:,1) )
444  axis(1:im,1:jm,1) = real_lony(ims:ime,jms:jme) / d2r
445  axis_name(1:2) = (/'x ','yh'/)
446  call historyputassociatedcoordinates( 'lon_xv', 'longitude (half level xv)', &
447  'degrees_east' , axis_name(1:2), axis(:,:,1) )
448  axis(1:im,1:jm,1) = real_lonxy(ims:ime,jms:jme) / d2r
449  axis_name(1:2) = (/'xh','yh'/)
450  call historyputassociatedcoordinates( 'lon_uv', 'longitude (half level uv)', &
451  'degrees_east' , axis_name(1:2), axis(:,:,1) )
452 
453  axis(1:im,1:jm,1) = real_lat(ims:ime,jms:jme) / d2r
454  axis_name(1:2) = (/'x ','y '/)
455  call historyputassociatedcoordinates( 'lat' , 'latitude' , &
456  'degrees_north', axis_name(1:2), axis(:,:,1) )
457 
458  axis(1:im,1:jm,1) = real_latx(ims:ime,jms:jme) / d2r
459  axis_name(1:2) = (/'xh','y '/)
460  call historyputassociatedcoordinates( 'lat_uy', 'latitude (half level uy)' , &
461  'degrees_north', axis_name(1:2), axis(:,:,1) )
462  axis(1:im,1:jm,1) = real_laty(ims:ime,jms:jme) / d2r
463  axis_name(1:2) = (/'x ','yh'/)
464  call historyputassociatedcoordinates( 'lat_xv', 'latitude (half level xv)' , &
465  'degrees_north', axis_name(1:2), axis(:,:,1) )
466  axis(1:im,1:jm,1) = real_latxy(ims:ime,jms:jme) / d2r
467  axis_name(1:2) = (/'xh','yh'/)
468  call historyputassociatedcoordinates( 'lat_uv', 'latitude (half level uv)' , &
469  'degrees_north', axis_name(1:2), axis(:,:,1) )
470 
471  axis(1:im,1:jm,1) = topo_zsfc(ims:ime,jms:jme)
472  axis_name(1:2) = (/'x ','y '/)
473  call historyputassociatedcoordinates( 'topo', 'topography', &
474  'm', axis_name(1:2), axis(:,:,1) )
475 
476  axis(1:im,1:jm,1) = landuse_frac_land(ims:ime,jms:jme)
477  axis_name(1:2) = (/'x ','y '/)
478  call historyputassociatedcoordinates( 'lsmask', 'fraction for land-sea mask', &
479  '0-1', axis_name(1:2), axis(:,:,1) )
480 
481  return
482  end subroutine hist_put_axes
483 
484  !-----------------------------------------------------------------------------
486  subroutine hist_reg( &
487  itemid, &
488  zinterp, &
489  item, &
490  desc, &
491  unit, &
492  ndim, &
493  xdim, &
494  ydim, &
495  zdim )
496  use gtool_history, only: &
498  use scale_time, only: &
499  nowstep => time_nowstep
500  implicit none
501 
502  integer, intent(out) :: itemid
503  logical, intent(out) :: zinterp
504  character(len=*), intent(in) :: item
505  character(len=*), intent(in) :: desc
506  character(len=*), intent(in) :: unit
507  integer, intent(in) :: ndim
508 
509  character(len=*), intent(in), optional :: xdim
510  character(len=*), intent(in), optional :: ydim
511  character(len=*), intent(in), optional :: zdim
512 
513  logical :: xh
514  logical :: yh
515 
516  logical :: existed
517 
518  character(len=H_SHORT) :: dims(3)
519  !---------------------------------------------------------------------------
520 
521  call prof_rapstart('FILE_O_NetCDF', 2)
522 
523  if ( ndim == 1 ) then
524 
525  dims(1) = "z"
526  if ( present(zdim) ) then
527  if ( zdim=='half' ) then
528  dims(1) = "zh"
529  end if
530  end if
531 
532  else
533 
534  xh = .false.
535  if ( present(xdim) ) then
536  if ( xdim=='half' ) xh = .true.
537  end if
538  yh = .false.
539  if ( present(ydim) ) then
540  if ( ydim=='half' ) yh = .true.
541  end if
542 
543  if ( xh .AND. yh ) then
544  dims(1) = 'lon_uv'
545  dims(2) = 'lat_uv'
546  dims(3) = 'height_uvz'
547  else if ( xh ) then
548  dims(1) = 'lon_uy'
549  dims(2) = 'lat_uy'
550  dims(3) = 'height_uyz'
551  else if ( yh ) then
552  dims(1) = 'lon_xv'
553  dims(2) = 'lat_xv'
554  dims(3) = 'height_xvz'
555  else
556  dims(1) = 'lon'
557  dims(2) = 'lat'
558  dims(3) = 'height'
559  end if
560 
561  if ( present(zdim) ) then
562  if ( zdim=='land' ) then
563  dims(3) = 'lz'
564  else if ( zdim=='landhalf' ) then
565  dims(3) = 'lzh'
566  else if ( zdim=='urban' ) then
567  dims(3) = 'uz'
568  else if ( zdim=='urbanhalf' ) then
569  dims(3) = 'uzh'
570  else if ( zdim=='half' ) then
571  if ( xh .AND. yh ) then
572  dims(3) = 'height_uvw'
573  else if ( xh ) then
574  dims(3) = 'height_uyw'
575  else if ( yh ) then
576  dims(3) = 'height_xvw'
577  else
578  dims(3) = 'height_xyw'
579  end if
580  end if
581  endif
582 
583  end if
584 
585  call historyaddvariable( item, & ! [IN]
586  dims(1:ndim), & ! [IN]
587  desc, & ! [IN]
588  unit, & ! [IN]
589  nowstep, & ! [IN]
590  itemid, & ! [OUT]
591  zinterp, & ! [OUT]
592  existed ) ! [OUT]
593 
594  call prof_rapend ('FILE_O_NetCDF', 2)
595 
596  return
597  end subroutine hist_reg
598 
599  !-----------------------------------------------------------------------------
601  subroutine hist_query( &
602  itemid, &
603  answer )
604  use gtool_history, only: &
606  use scale_time, only: &
608  implicit none
609 
610  integer, intent(in) :: itemid
611  logical, intent(out) :: answer
612 
613  !---------------------------------------------------------------------------
614 
615  answer = .false.
616 
617  if ( .NOT. enabled ) return
618 
619  if ( itemid < 0 ) return
620 
621  call prof_rapstart('FILE_O_NetCDF', 2)
622 
623  call historyquery(itemid, time_nowstep, answer)
624 
625  call prof_rapend ('FILE_O_NetCDF', 2)
626 
627  return
628  end subroutine hist_query
629 
630  !-----------------------------------------------------------------------------
632  subroutine hist_put_0d( &
633  itemid, &
634  var )
635  use gtool_history, only: &
636  historyput
637  use scale_time, only: &
639  implicit none
640 
641  integer, intent(in) :: itemid
642  real(RP), intent(in) :: var
643 
644  !---------------------------------------------------------------------------
645 
646  if ( .NOT. enabled ) return
647 
648  if ( itemid < 0 ) return
649 
650  call prof_rapstart('FILE_O_NetCDF', 2)
651 
652  call historyput(itemid, time_nowstep, var)
653 
654  call prof_rapend ('FILE_O_NetCDF', 2)
655 
656  return
657  end subroutine hist_put_0d
658 
659  !-----------------------------------------------------------------------------
661  subroutine hist_put_1d( &
662  itemid, &
663  var )
664  use gtool_history, only: &
665  historyput
666  use scale_time, only: &
668  implicit none
669 
670  integer, intent(in) :: itemid
671  real(RP), intent(in) :: var(:)
672 
673  real(RP) :: var2(kmax)
674  integer :: k
675  !---------------------------------------------------------------------------
676 
677  if ( .NOT. enabled ) return
678 
679  if ( itemid < 0 ) return
680 
681  call prof_rapstart('FILE_O_NetCDF', 2)
682 
683  do k = 1, kmax
684  var2(k) = var(ks+k-1)
685  enddo
686 
687  call historyput(itemid, time_nowstep, var2)
688 
689  call prof_rapend ('FILE_O_NetCDF', 2)
690 
691  return
692  end subroutine hist_put_1d
693 
694  !-----------------------------------------------------------------------------
696  subroutine hist_put_2d( &
697  itemid, &
698  var, &
699  nohalo )
700  use gtool_file, only: &
701  rmiss
702  use gtool_history, only: &
703  historyput
704  use scale_time, only: &
706  implicit none
707 
708  integer, intent(in) :: itemid
709  real(RP), intent(in) :: var(:,:)
710 
711  logical, intent(in), optional :: nohalo
712 
713  real(RP) :: var2(im*jm)
714  integer :: i, j
715  logical :: nohalo_
716  !---------------------------------------------------------------------------
717 
718  if ( .NOT. enabled ) return
719 
720  if ( itemid < 0 ) return
721 
722  call prof_rapstart('FILE_O_NetCDF', 2)
723 
724  do j = 1, jm
725  do i = 1, im
726  var2(i + (j-1)*im) = var(ims+i-1,jms+j-1)
727  enddo
728  enddo
729 
730  nohalo_ = .false.
731  if ( present(nohalo) ) nohalo_ = nohalo
732 
733  if ( nohalo_ ) then
734  ! W halo
735  do j = 1, jm
736  do i = 1, is-ims
737  var2(i + (j-1)*im) = rmiss
738  enddo
739  enddo
740  ! E halo
741  do j = 1, jm
742  do i = ie-ims+2, ime-ims+1
743  var2(i + (j-1)*im) = rmiss
744  enddo
745  enddo
746  ! S halo
747  do j = 1, js-jms
748  do i = 1, im
749  var2(i + (j-1)*im) = rmiss
750  enddo
751  enddo
752  ! N halo
753  do j = je-jms+2, jme-jms+1
754  do i = 1, im
755  var2(i + (j-1)*im) = rmiss
756  enddo
757  enddo
758  end if
759 
760  call historyput(itemid, time_nowstep, var2)
761 
762  call prof_rapend ('FILE_O_NetCDF', 2)
763 
764  return
765  end subroutine hist_put_2d
766 
767  !-----------------------------------------------------------------------------
769  subroutine hist_put_3d( &
770  itemid, &
771  var, &
772  zinterp, &
773  xdim, &
774  ydim, &
775  zdim, &
776  nohalo )
777  use gtool_file, only: &
778  rmiss
779  use gtool_history, only: &
780  historyput
781  use scale_time, only: &
783  use scale_interpolation, only: &
786  implicit none
787 
788  integer, intent(in) :: itemid
789  real(RP), intent(in) :: var(:,:,:)
790  logical, intent(in) :: zinterp
791 
792  character(len=*), intent(in), optional :: xdim
793  character(len=*), intent(in), optional :: ydim
794  character(len=*), intent(in), optional :: zdim
795  logical, intent(in), optional :: nohalo
796 
797  intrinsic shape
798  integer :: s(3)
799 
800  character(len=H_SHORT) :: xd, yd, zd
801 
802  real(RP), allocatable :: var_Z(:,:,:)
803  real(RP), allocatable :: var2 (:)
804 
805  integer :: i, j, k
806  integer :: isize, jsize, ksize
807  integer :: iall, jall, kall
808  integer :: istart, jstart, kstart
809 
810  logical :: nohalo_
811  !---------------------------------------------------------------------------
812 
813  if ( .NOT. enabled ) return
814 
815  if ( itemid < 0 ) return
816 
817  call prof_rapstart('FILE_O_NetCDF', 2)
818 
819  xd = ''
820  yd = ''
821  zd = ''
822  if( present(xdim) ) xd = xdim
823  if( present(ydim) ) yd = ydim
824  if( present(zdim) ) zd = zdim
825 
826  nohalo_ = .false.
827  if ( present(nohalo) ) nohalo_ = nohalo
828 
829 
830  ! select dimension
831  select case ( xd )
832  case default
833  isize = im
834  iall = ia
835  istart = ims
836  end select
837  select case ( yd )
838  case default
839  jsize = jm
840  jall = ja
841  jstart = jms
842  end select
843  select case ( zd )
844  case ('land')
845  ksize = lkmax
846  kall = lke
847  kstart = lks
848  case ('urban')
849  ksize = ukmax
850  kall = uke
851  kstart = uks
852  case default
853  ksize = kmax
854  kall = ka
855  kstart = ks
856  end select
857 
858  allocate( var_z(kall,iall,jall) )
859  allocate( var2(ksize*isize*jsize) )
860 
861  s = shape(var)
862  if ( s(1) == 1 ) then
863 
864  do j = 1, jsize
865  do i = 1, isize
866  var2(i + (j-1)*isize) = var(1,istart+i-1,jstart+j-1)
867  enddo
868  enddo
869 
870  if ( nohalo_ ) then
871  ! W halo
872  do j = 1, jsize
873  do i = 1, is-ims
874  var2(i + (j-1)*isize) = rmiss
875  enddo
876  enddo
877  ! E halo
878  do j = 1, jsize
879  do i = ie-ims+2, ime-ims+1
880  var2(i + (j-1)*isize) = rmiss
881  enddo
882  enddo
883  ! S halo
884  do j = 1, js-jms
885  do i = 1, isize
886  var2(i + (j-1)*isize) = rmiss
887  enddo
888  enddo
889  ! N halo
890  do j = je-jms+2, jme-jms+1
891  do i = 1, isize
892  var2(i + (j-1)*isize) = rmiss
893  enddo
894  enddo
895  end if
896 
897  call historyput(itemid, time_nowstep, var2(1:isize*jsize))
898 
899  else
900  if ( ksize == kmax &
901  .AND. zinterp &
902  .AND. interp_available ) then
903  call prof_rapstart('FILE_O_interp', 2)
904  call interp_vertical_xi2z( var(:,:,:), & ! [IN]
905  var_z(:,:,:) ) ! [OUT]
906  call prof_rapend ('FILE_O_interp', 2)
907 
908  do k = 1, ksize
909  do j = 1, jsize
910  do i = 1, isize
911  var2(i + (j-1)*isize + (k-1)*jsize*isize) = var_z(kstart+k-1,istart+i-1,jstart+j-1)
912  enddo
913  enddo
914  enddo
915  else
916  do k = 1, ksize
917  do j = 1, jsize
918  do i = 1, isize
919  var2(i + (j-1)*isize + (k-1)*jsize*isize) = var(kstart+k-1,istart+i-1,jstart+j-1)
920  enddo
921  enddo
922  enddo
923  endif
924 
925  if ( nohalo_ ) then
926  ! W halo
927  do k = 1, ksize
928  do j = 1, jsize
929  do i = 1, is-ims
930  var2(i + (j-1)*isize + (k-1)*jsize*isize) = rmiss
931  enddo
932  enddo
933  enddo
934  ! E halo
935  do k = 1, ksize
936  do j = 1, jsize
937  do i = ie-ims+2, ime-ims+1
938  var2(i + (j-1)*isize + (k-1)*jsize*isize) = rmiss
939  enddo
940  enddo
941  enddo
942  ! S halo
943  do k = 1, ksize
944  do j = 1, js-jms
945  do i = 1, isize
946  var2(i + (j-1)*isize + (k-1)*jsize*isize) = rmiss
947  enddo
948  enddo
949  enddo
950  ! N halo
951  do k = 1, ksize
952  do j = je-jms+2, jme-jms+1
953  do i = 1, isize
954  var2(i + (j-1)*isize + (k-1)*jsize*isize) = rmiss
955  enddo
956  enddo
957  enddo
958  end if
959 
960  call historyput(itemid, time_nowstep, var2)
961 
962  endif
963 
964  deallocate( var_z )
965  deallocate( var2 )
966 
967  call prof_rapend ('FILE_O_NetCDF', 2)
968 
969  return
970  end subroutine hist_put_3d
971 
972  !-----------------------------------------------------------------------------
974  subroutine hist_in_0d( &
975  var, &
976  item, &
977  desc, &
978  unit )
979  implicit none
980 
981  real(RP), intent(in) :: var
982  character(len=*), intent(in) :: item
983  character(len=*), intent(in) :: desc
984  character(len=*), intent(in) :: unit
985 
986  integer :: itemid
987  logical :: zinterp
988  logical :: do_put
989  !---------------------------------------------------------------------------
990 
991  if ( .NOT. enabled ) return
992 
993 
994  call hist_reg( itemid, & ! [OUT]
995  zinterp, & ! [OUT]
996  item, desc, unit, 0 ) ! [IN]
997 
998 
999  call hist_query( itemid, & ! [IN]
1000  do_put ) ! [OUT]
1001 
1002  if ( do_put ) then
1003  call hist_put( itemid, var ) ! [IN]
1004  endif
1005 
1006  return
1007  end subroutine hist_in_0d
1008 
1009  !-----------------------------------------------------------------------------
1011  subroutine hist_in_1d( &
1012  var, &
1013  item, &
1014  desc, &
1015  unit, &
1016  zdim )
1017  implicit none
1018 
1019  real(RP), intent(in) :: var(:)
1020  character(len=*), intent(in) :: item
1021  character(len=*), intent(in) :: desc
1022  character(len=*), intent(in) :: unit
1023 
1024  character(len=*), intent(in), optional :: zdim
1025 
1026  character(len=H_SHORT) :: zd
1027  integer :: itemid
1028  logical :: zinterp
1029  logical :: do_put
1030  !---------------------------------------------------------------------------
1031 
1032  if ( .NOT. enabled ) return
1033 
1034  zd = ''
1035  if( present(zdim) ) zd = zdim
1036 
1037  call hist_reg( itemid, & ! [OUT]
1038  zinterp, & ! [OUT]
1039  item, desc, unit, 1, & ! [IN]
1040  zdim = zd ) ! [IN]
1041 
1042  call hist_query( itemid, & ! [IN]
1043  do_put ) ! [OUT]
1044 
1045  if ( do_put ) then
1046  call hist_put( itemid, var ) ! [IN]
1047  endif
1048 
1049  return
1050  end subroutine hist_in_1d
1051 
1052  !-----------------------------------------------------------------------------
1054  subroutine hist_in_2d( &
1055  var, &
1056  item, &
1057  desc, &
1058  unit, &
1059  xdim, &
1060  ydim, &
1061  nohalo )
1062  implicit none
1063 
1064  real(RP), intent(in) :: var(:,:)
1065  character(len=*), intent(in) :: item
1066  character(len=*), intent(in) :: desc
1067  character(len=*), intent(in) :: unit
1068 
1069  character(len=*), intent(in), optional :: xdim
1070  character(len=*), intent(in), optional :: ydim
1071  logical, intent(in), optional :: nohalo
1072 
1073  character(len=H_SHORT) :: xd, yd
1074  integer :: itemid
1075  logical :: zinterp
1076  logical :: do_put
1077  !---------------------------------------------------------------------------
1078 
1079  if ( .NOT. enabled ) return
1080 
1081  xd = ''
1082  yd = ''
1083  if( present(xdim) ) xd = xdim
1084  if( present(ydim) ) yd = ydim
1085 
1086  call hist_reg( itemid, & ! [OUT]
1087  zinterp, & ! [OUT]
1088  item, desc, unit, 2, & ! [IN]
1089  xdim = xd, ydim = yd ) ! [IN]
1090 
1091  call hist_query( itemid, & ! [IN]
1092  do_put ) ! [OUT]
1093 
1094  if ( do_put ) then
1095  call hist_put( itemid, var, nohalo = nohalo ) ! [IN]
1096  endif
1097 
1098  return
1099  end subroutine hist_in_2d
1100 
1101  !-----------------------------------------------------------------------------
1103  subroutine hist_in_3d( &
1104  var, &
1105  item, &
1106  desc, &
1107  unit, &
1108  xdim, &
1109  ydim, &
1110  zdim, &
1111  nohalo )
1112  implicit none
1113 
1114  real(RP), intent(in) :: var(:,:,:)
1115  character(len=*), intent(in) :: item
1116  character(len=*), intent(in) :: desc
1117  character(len=*), intent(in) :: unit
1118 
1119  character(len=*), intent(in), optional :: xdim
1120  character(len=*), intent(in), optional :: ydim
1121  character(len=*), intent(in), optional :: zdim
1122  logical, intent(in), optional :: nohalo
1123 
1124  character(len=H_SHORT) :: xd, yd, zd
1125  integer :: itemid
1126  logical :: zinterp
1127  logical :: do_put
1128  !---------------------------------------------------------------------------
1129 
1130  if ( .NOT. enabled ) return
1131 
1132  xd = ''
1133  yd = ''
1134  zd = ''
1135  if( present(xdim) ) xd = xdim
1136  if( present(ydim) ) yd = ydim
1137  if( present(zdim) ) zd = zdim
1138 
1139  call hist_reg( itemid, & ! [OUT]
1140  zinterp, & ! [OUT]
1141  item, desc, unit, 3, & ! [IN]
1142  xdim = xd, ydim = yd, zdim=zd ) ! [IN]
1143 
1144  call hist_query( itemid, & ! [IN]
1145  do_put ) ! [OUT]
1146 
1147  if ( do_put ) then
1148  call hist_put( itemid, var, zinterp, & ! [IN]
1149  xdim = xd, ydim = yd, zdim=zd, & ! [IN]
1150  nohalo = nohalo ) ! [IN]
1151  endif
1152 
1153  return
1154  end subroutine hist_in_3d
1155 
1156  !-----------------------------------------------------------------------------
1158  subroutine hist_get_1d( &
1159  var, &
1160  basename, &
1161  varname, &
1162  step, &
1163  allow_missing )
1164  use gtool_history, only: &
1165  historyget
1166  implicit none
1167 
1168  real(RP), intent(out) :: var(:)
1169  character(len=*), intent(in) :: basename
1170  character(len=*), intent(in) :: varname
1171  integer, intent(in) :: step
1172 
1173  logical, intent(in), optional :: allow_missing
1174 
1175  logical :: am
1176  !---------------------------------------------------------------------------
1177 
1178  call prof_rapstart('FILE_I_NetCDF', 2)
1179 
1180  am = .false.
1181  if( present(allow_missing) ) am = allow_missing
1182 
1183  call historyget( var(:), & ! [OUT]
1184  basename, & ! [IN]
1185  varname, & ! [IN]
1186  step, & ! [IN]
1187  allow_missing=am ) ! [IN]
1188 
1189  call prof_rapend ('FILE_I_NetCDF', 2)
1190 
1191  return
1192  end subroutine hist_get_1d
1193 
1194  !-----------------------------------------------------------------------------
1196  subroutine hist_get_2d( &
1197  var, &
1198  basename, &
1199  varname, &
1200  step, &
1201  allow_missing )
1202  use gtool_history, only: &
1203  historyget
1204  implicit none
1205 
1206  real(RP), intent(out) :: var(:,:)
1207  character(len=*), intent(in) :: basename
1208  character(len=*), intent(in) :: varname
1209  integer, intent(in) :: step
1210 
1211  logical, intent(in), optional :: allow_missing
1212 
1213  logical :: am
1214  !---------------------------------------------------------------------------
1215 
1216  call prof_rapstart('FILE_I_NetCDF', 2)
1217 
1218  am = .false.
1219  if( present(allow_missing) ) am = allow_missing
1220 
1221  call historyget( var(:,:), & ! [OUT]
1222  basename, & ! [IN]
1223  varname, & ! [IN]
1224  step, & ! [IN]
1225  allow_missing=am ) ! [IN]
1226 
1227  call prof_rapend ('FILE_I_NetCDF', 2)
1228 
1229  return
1230  end subroutine hist_get_2d
1231 
1232  !-----------------------------------------------------------------------------
1234  subroutine hist_get_3d( &
1235  var, &
1236  basename, &
1237  varname, &
1238  step, &
1239  allow_missing )
1240  use gtool_history, only: &
1241  historyget
1242  implicit none
1243 
1244  real(RP), intent(out) :: var(:,:,:)
1245  character(len=*), intent(in) :: basename
1246  character(len=*), intent(in) :: varname
1247  integer, intent(in) :: step
1248 
1249  logical, intent(in), optional :: allow_missing
1250 
1251  logical :: am
1252  !---------------------------------------------------------------------------
1253 
1254  call prof_rapstart('FILE_I_NetCDF', 2)
1255 
1256  am = .false.
1257  if( present(allow_missing) ) am = allow_missing
1258 
1259  call historyget( var(:,:,:), & ! [OUT]
1260  basename, & ! [IN]
1261  varname, & ! [IN]
1262  step, & ! [IN]
1263  allow_missing=am ) ! [IN]
1264 
1265  call prof_rapend ('FILE_I_NetCDF', 2)
1266 
1267  return
1268  end subroutine hist_get_3d
1269 
1270  !-----------------------------------------------------------------------------
1272  subroutine hist_write
1273  use gtool_history, only: &
1275  use scale_time, only: &
1276  time_nowstep
1277  implicit none
1278  !---------------------------------------------------------------------------
1279 
1280  call prof_rapstart('FILE_O_NetCDF', 2)
1281 
1282  call historywriteall( time_nowstep ) ![IN]
1283 
1284  call prof_rapend ('FILE_O_NetCDF', 2)
1285 
1286  return
1287  end subroutine hist_write
1288 
1289 end module scale_history
1290 !-------------------------------------------------------------------------------
integer, public imax
of computational cells: x
integer, public time_nowstep
current step [number]
Definition: scale_time.F90:70
integer, public is
start point of inner domain: x, local
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public je
end point of inner domain: y, local
subroutine, public hist_reg(itemid, zinterp, item, desc, unit, ndim, xdim, ydim, zdim)
Register/Append variable to history file.
real(rp), dimension(:), allocatable, public grid_cyg
center coordinate [m]: y, global
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:), allocatable, public grid_cbfyg
center buffer factor [0-1]: y, global
integer, public jeb
real(dp), public time_nowms
subsecond part of current time [millisec]
Definition: scale_time.F90:66
real(dp), parameter, public rmiss
Definition: gtool_file.f90:110
real(rp), dimension(:), allocatable, public grid_cxg
center coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fdy
y-length of grid(j+1) to grid(j) [m]
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public hist_switch(switch)
set switch
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fxg
face coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fbfy
face buffer factor [0-1]: y
module STDIO
Definition: scale_stdio.F90:12
real(dp), public time_startdaysec
second of start time [sec]
Definition: scale_time.F90:74
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
integer, public imaxb
subroutine, public interp_vertical_xi2z(var, var_Z)
Reset random seed.
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
module GRID (cartesian) for land
subroutine hist_put_1d(itemid, var)
Put 1D data to history buffer.
subroutine hist_put_3d(itemid, var, zinterp, xdim, ydim, zdim, nohalo)
Put 3D data to history buffer.
integer, public jmaxb
module Gtool_History
real(rp), dimension(:), allocatable, public grid_fx
face coordinate [m]: x, local
real(rp), dimension(:), allocatable, public grid_lfz
face coordinate [m]: z, local=global
integer, public ieb
subroutine hist_in_0d(var, item, desc, unit)
Wrapper routine of HIST_reg+HIST_put 0D.
module grid index
real(rp), dimension(:), allocatable, public grid_ucdz
z-length of control volume [m]
subroutine, public hist_query(itemid, answer)
Check time to putting data.
integer, public ia
of x whole cells (local, with HALO)
character(len=h_mid), public h_source
for file header
Definition: scale_stdio.F90:50
real(rp), dimension(:,:), allocatable, public real_latx
latitude at staggered point (uy) [rad,-pi,pi]
module GRID (cartesian) for urban
real(rp), dimension(:), allocatable, public grid_fbfx
face buffer factor [0-1]: x
real(rp), dimension(:), allocatable, public grid_fdz
z-length of grid(k+1) to grid(k) [m]
module GRID (real space)
module LANDUSE
real(rp), dimension(:), allocatable, public grid_ufz
face coordinate [m]: z, local=global
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:36
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fbfxg
face buffer factor [0-1]: x, global
subroutine, public hist_write
Flush history buffer to file.
integer, public time_offset_year
time offset [year]
Definition: scale_time.F90:73
integer, public kmax
of computational cells: z
real(rp), dimension(:), allocatable, public grid_fz
face coordinate [m]: z, local=global
real(rp), dimension(:,:), allocatable, public real_lonxy
longitude at staggered point (uv) [rad,0-2pi]
subroutine hist_get_3d(var, basename, varname, step, allow_missing)
Get 3D data from file.
real(rp), dimension(:), allocatable, public grid_fbfz
face buffer factor [0-1]: z
subroutine hist_get_2d(var, basename, varname, step, allow_missing)
Get 2D data from file.
subroutine, public historyaddvariable(varname, dims, desc, units, now_step, id, zinterp, existed, options)
subroutine, public historyinit(title, source, institution, array_size, master, myrank, rankidx, time_start, time_interval, time_units, time_since, default_basename, default_tinterval, default_tunit, default_taverage, default_zinterp, default_datatype, namelist_filename, namelist_fid)
integer, public js
start point of inner domain: y, local
real(rp), dimension(:), allocatable, public grid_cbfx
center buffer factor [0-1]: x
module TIME
Definition: scale_time.F90:15
module PROCESS
real(rp), dimension(:), allocatable, public grid_fbfyg
face buffer factor [0-1]: y, global
real(rp), dimension(:), allocatable, public grid_ucz
center coordinate [m]: z, local=global
module CONSTANT
Definition: scale_const.F90:14
integer, parameter, public prc_masterrank
master process in each communicator
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:), allocatable, public grid_cbfz
center buffer factor [0-1]: z
integer, public prc_myrank
process num in local communicator
real(rp), dimension(:), allocatable, public grid_cx
center coordinate [m]: x, local
module GRID (cartesian)
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
module RM PROCESS
subroutine, public hist_setup
Setup.
subroutine hist_put_0d(itemid, var)
Put 1D data to history buffer.
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:), allocatable, public grid_lcdz
z-length of control volume [m]
character(len=h_mid), public h_institute
for file header
Definition: scale_stdio.F90:51
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:,:), allocatable, public real_lon
longitude [rad,0-2pi]
logical, public interp_available
vertical interpolation is available? (have meaning?)
module PRECISION
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
subroutine, public historywriteall(step_now)
real(rp), dimension(:), allocatable, public grid_fdx
x-length of grid(i+1) to grid(i) [m]
real(rp), dimension(:), allocatable, public grid_lcz
center coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_cdy
y-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_cbfy
center buffer factor [0-1]: y
module HISTORY
module TOPOGRAPHY
integer, public isb
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
module land grid index
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:,:), allocatable, public real_lat
latitude [rad,-pi,pi]
real(rp), dimension(:,:), allocatable, public real_lony
longitude at staggered point (xv) [rad,0-2pi]
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public jsb
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
real(rp), dimension(:,:), allocatable, public real_latxy
latitude at staggered point (uv) [rad,-pi,pi]
real(rp), dimension(:,:), allocatable, public real_lonx
longitude at staggered point (uy) [rad,0-2pi]
real(rp), dimension(:), allocatable, public grid_cdx
x-length of control volume [m]
real(rp), dimension(:,:), allocatable, public landuse_frac_land
land fraction
integer, public jmax
of computational cells: y
subroutine, public historyquery(itemid, step_next, answer)
real(rp), dimension(:), allocatable, public grid_fyg
face coordinate [m]: y, global
real(rp), dimension(:), allocatable, public grid_cy
center coordinate [m]: y, local
subroutine hist_put_2d(itemid, var, nohalo)
Put 2D data to history buffer.
module INTERPOLATION
real(rp), dimension(:), allocatable, public grid_cbfxg
center buffer factor [0-1]: x, global
module urban grid index
real(rp), dimension(:,:), allocatable, public real_laty
latitude at staggered point (xv) [rad,-pi,pi]
integer, public ja
of y whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local