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