45 module procedure hist_in_1d
46 module procedure hist_in_2d
47 module procedure hist_in_3d
54 end interface hist_put
56 module procedure hist_get_1d
59 end interface hist_get
69 private :: hist_put_axes
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 =
'' 104 logical :: HIST_BND = .true.
106 namelist / param_hist / &
109 real(DP) :: start_daysec
110 integer :: rankidx(2)
115 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[HISTORY] / Categ[ATMOS-RM IO] / Origin[SCALElib]' 121 if(
io_l )
write(
io_fid_log,*)
'*** Not found namelist. Default used.' 122 elseif( ierr > 0 )
then 123 write(*,*)
'xxx Not appropriate names in namelist PARAM_HIST. Check!' 140 write(history_t_since,
'(i4.4,a1,i2.2,a1,i2.2,a1,i2.2,a1,i2.2,a1,i2.2)') &
171 time_units = history_t_units, &
172 time_since = history_t_since, &
187 logical,
intent(in) :: switch
196 subroutine hist_put_axes
199 historyputassociatedcoordinates
254 real(RP) :: AXIS (im,jm,
kmax)
255 character(len=2) :: AXIS_name(3)
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) )
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. )
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. )
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 )
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 )
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 )
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 )
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)
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)
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)
314 axis(1:im,1:jm,k) =
real_cz(k+
ks-1,ims:ime,jms:jme)
316 axis_name(1:3) = (/
'x ',
'y ',
'z '/)
317 call historyputassociatedcoordinates(
'height' ,
'height' , &
318 'm' , axis_name(1:3), axis(:,:,:) )
321 axis(1:im,1:jm,k) =
real_fz(k+
ks-1,ims:ime,jms:jme)
323 axis_name(1:3) = (/
'x ',
'y ',
'zh'/)
324 call historyputassociatedcoordinates(
'height_xyw' ,
'height (half level xyw)' , &
325 'm' , axis_name(1:3), axis(:,:,:) )
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
336 axis(im,j,k) =
real_cz(k+
ks-1,ims+im-1,jms+j-1)
339 axis_name(1:3) = (/
'xh',
'y ',
'z '/)
340 call historyputassociatedcoordinates(
'height_uyz' ,
'height (half level uyz)' , &
341 'm' , axis_name(1:3), axis(:,:,:) )
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
352 axis(i,jm,k) =
real_cz(k+
ks-1,ims+i-1,jms+jm-1)
355 axis_name(1:3) = (/
'x ',
'yh',
'z '/)
356 call historyputassociatedcoordinates(
'height_xvz' ,
'height (half level xvz)' , &
357 'm' , axis_name(1:3), axis(:,:,:) )
369 axis(i,jm,k) =
real_cz(k+
ks-1,ims+i-1,jms+jm-1)
374 axis(im,j,k) =
real_cz(k+
ks-1,ims+im-1,jms+j-1)
377 axis_name(1:3) = (/
'xh',
'yh',
'z '/)
378 call historyputassociatedcoordinates(
'height_uvz' ,
'height (half level uvz)' , &
379 'm' , axis_name(1:3), axis(:,:,:) )
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
390 axis(im,j,k) =
real_fz(k+
ks-1,ims+im-1,jms+j-1)
393 axis_name(1:3) = (/
'xh',
'y ',
'zh'/)
394 call historyputassociatedcoordinates(
'height_uyw' ,
'height (half level uyw)' , &
395 'm' , axis_name(1:3), axis(:,:,:) )
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
406 axis(i,jm,k) =
real_fz(k+
ks-1,ims+i-1,jms+jm-1)
409 axis_name(1:3) = (/
'x ',
'yh',
'zh'/)
410 call historyputassociatedcoordinates(
'height_xvw' ,
'height (half level xvw)' , &
411 'm' , axis_name(1:3), axis(:,:,:) )
423 axis(i,jm,k) =
real_fz(k+
ks-1,ims+i-1,jms+jm-1)
428 axis(im,j,k) =
real_fz(k+
ks-1,ims+im-1,jms+j-1)
431 axis_name(1:3) = (/
'xh',
'yh',
'zh'/)
432 call historyputassociatedcoordinates(
'height_uvw' ,
'height (half level uvw)' , &
433 'm' , axis_name(1:3), axis(:,:,:) )
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) )
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) )
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) )
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) )
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) )
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) )
482 end subroutine hist_put_axes
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
509 character(len=*),
intent(in),
optional :: xdim
510 character(len=*),
intent(in),
optional :: ydim
511 character(len=*),
intent(in),
optional :: zdim
518 character(len=H_SHORT) :: dims(3)
523 if ( ndim == 1 )
then 526 if (
present(zdim) )
then 527 if ( zdim==
'half' )
then 535 if (
present(xdim) )
then 536 if ( xdim==
'half' ) xh = .true.
539 if (
present(ydim) )
then 540 if ( ydim==
'half' ) yh = .true.
543 if ( xh .AND. yh )
then 546 dims(3) =
'height_uvz' 550 dims(3) =
'height_uyz' 554 dims(3) =
'height_xvz' 561 if (
present(zdim) )
then 562 if ( zdim==
'land' )
then 564 else if ( zdim==
'landhalf' )
then 566 else if ( zdim==
'urban' )
then 568 else if ( zdim==
'urbanhalf' )
then 570 else if ( zdim==
'half' )
then 571 if ( xh .AND. yh )
then 572 dims(3) =
'height_uvw' 574 dims(3) =
'height_uyw' 576 dims(3) =
'height_xvw' 578 dims(3) =
'height_xyw' 610 integer,
intent(in) :: itemid
611 logical,
intent(out) :: answer
617 if ( .NOT. enabled )
return 619 if ( itemid < 0 )
return 641 integer,
intent(in) :: itemid
642 real(RP),
intent(in) :: var
646 if ( .NOT. enabled )
return 648 if ( itemid < 0 )
return 670 integer,
intent(in) :: itemid
671 real(RP),
intent(in) :: var(:)
673 real(RP) :: var2(
kmax)
677 if ( .NOT. enabled )
return 679 if ( itemid < 0 )
return 684 var2(k) = var(
ks+k-1)
708 integer,
intent(in) :: itemid
709 real(RP),
intent(in) :: var(:,:)
711 logical,
intent(in),
optional :: nohalo
713 real(RP) :: var2(im*jm)
718 if ( .NOT. enabled )
return 720 if ( itemid < 0 )
return 726 var2(i + (j-1)*im) = var(ims+i-1,jms+j-1)
731 if (
present(nohalo) ) nohalo_ = nohalo
737 var2(i + (j-1)*im) =
rmiss 742 do i =
ie-ims+2, ime-ims+1
743 var2(i + (j-1)*im) =
rmiss 749 var2(i + (j-1)*im) =
rmiss 753 do j =
je-jms+2, jme-jms+1
755 var2(i + (j-1)*im) =
rmiss 788 integer,
intent(in) :: itemid
789 real(RP),
intent(in) :: var(:,:,:)
790 logical,
intent(in) :: zinterp
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
800 character(len=H_SHORT) :: xd, yd, zd
802 real(RP),
allocatable :: var_Z(:,:,:)
803 real(RP),
allocatable :: var2 (:)
806 integer :: isize, jsize, ksize
807 integer :: iall, jall, kall
808 integer :: istart, jstart, kstart
813 if ( .NOT. enabled )
return 815 if ( itemid < 0 )
return 822 if(
present(xdim) ) xd = xdim
823 if(
present(ydim) ) yd = ydim
824 if(
present(zdim) ) zd = zdim
827 if (
present(nohalo) ) nohalo_ = nohalo
858 allocate( var_z(kall,iall,jall) )
859 allocate( var2(ksize*isize*jsize) )
862 if ( s(1) == 1 )
then 866 var2(i + (j-1)*isize) = var(1,istart+i-1,jstart+j-1)
874 var2(i + (j-1)*isize) =
rmiss 879 do i =
ie-ims+2, ime-ims+1
880 var2(i + (j-1)*isize) =
rmiss 886 var2(i + (j-1)*isize) =
rmiss 890 do j =
je-jms+2, jme-jms+1
892 var2(i + (j-1)*isize) =
rmiss 897 call historyput(itemid,
time_nowstep, var2(1:isize*jsize))
911 var2(i + (j-1)*isize + (k-1)*jsize*isize) = var_z(kstart+k-1,istart+i-1,jstart+j-1)
919 var2(i + (j-1)*isize + (k-1)*jsize*isize) = var(kstart+k-1,istart+i-1,jstart+j-1)
930 var2(i + (j-1)*isize + (k-1)*jsize*isize) =
rmiss 937 do i =
ie-ims+2, ime-ims+1
938 var2(i + (j-1)*isize + (k-1)*jsize*isize) =
rmiss 946 var2(i + (j-1)*isize + (k-1)*jsize*isize) =
rmiss 952 do j =
je-jms+2, jme-jms+1
954 var2(i + (j-1)*isize + (k-1)*jsize*isize) =
rmiss 981 real(RP),
intent(in) :: var
982 character(len=*),
intent(in) :: item
983 character(len=*),
intent(in) :: desc
984 character(len=*),
intent(in) :: unit
991 if ( .NOT. enabled )
return 996 item, desc, unit, 0 )
1003 call hist_put( itemid, var )
1011 subroutine hist_in_1d( &
1019 real(RP),
intent(in) :: var(:)
1020 character(len=*),
intent(in) :: item
1021 character(len=*),
intent(in) :: desc
1022 character(len=*),
intent(in) :: unit
1024 character(len=*),
intent(in),
optional :: zdim
1026 character(len=H_SHORT) :: zd
1032 if ( .NOT. enabled )
return 1035 if(
present(zdim) ) zd = zdim
1039 item, desc, unit, 1, &
1046 call hist_put( itemid, var )
1050 end subroutine hist_in_1d
1054 subroutine hist_in_2d( &
1064 real(RP),
intent(in) :: var(:,:)
1065 character(len=*),
intent(in) :: item
1066 character(len=*),
intent(in) :: desc
1067 character(len=*),
intent(in) :: unit
1069 character(len=*),
intent(in),
optional :: xdim
1070 character(len=*),
intent(in),
optional :: ydim
1071 logical,
intent(in),
optional :: nohalo
1073 character(len=H_SHORT) :: xd, yd
1079 if ( .NOT. enabled )
return 1083 if(
present(xdim) ) xd = xdim
1084 if(
present(ydim) ) yd = ydim
1088 item, desc, unit, 2, &
1089 xdim = xd, ydim = yd )
1095 call hist_put( itemid, var, nohalo = nohalo )
1099 end subroutine hist_in_2d
1103 subroutine hist_in_3d( &
1114 real(RP),
intent(in) :: var(:,:,:)
1115 character(len=*),
intent(in) :: item
1116 character(len=*),
intent(in) :: desc
1117 character(len=*),
intent(in) :: unit
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
1124 character(len=H_SHORT) :: xd, yd, zd
1130 if ( .NOT. enabled )
return 1135 if(
present(xdim) ) xd = xdim
1136 if(
present(ydim) ) yd = ydim
1137 if(
present(zdim) ) zd = zdim
1141 item, desc, unit, 3, &
1142 xdim = xd, ydim = yd, zdim=zd )
1148 call hist_put( itemid, var, zinterp, &
1149 xdim = xd, ydim = yd, zdim=zd, &
1154 end subroutine hist_in_3d
1158 subroutine hist_get_1d( &
1168 real(RP),
intent(out) :: var(:)
1169 character(len=*),
intent(in) :: basename
1170 character(len=*),
intent(in) :: varname
1171 integer,
intent(in) :: step
1173 logical,
intent(in),
optional :: allow_missing
1181 if(
present(allow_missing) ) am = allow_missing
1183 call historyget( var(:), &
1192 end subroutine hist_get_1d
1206 real(RP),
intent(out) :: var(:,:)
1207 character(len=*),
intent(in) :: basename
1208 character(len=*),
intent(in) :: varname
1209 integer,
intent(in) :: step
1211 logical,
intent(in),
optional :: allow_missing
1219 if(
present(allow_missing) ) am = allow_missing
1221 call historyget( var(:,:), &
1244 real(RP),
intent(out) :: var(:,:,:)
1245 character(len=*),
intent(in) :: basename
1246 character(len=*),
intent(in) :: varname
1247 integer,
intent(in) :: step
1249 logical,
intent(in),
optional :: allow_missing
1257 if(
present(allow_missing) ) am = allow_missing
1259 call historyget( var(:,:,:), &
integer, public imax
of computational cells: x
integer, public time_nowstep
current step [number]
integer, public is
start point of inner domain: x, local
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
real(dp), public time_nowms
subsecond part of current time [millisec]
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)
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
real(dp), public time_startdaysec
second of start time [sec]
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
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
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.
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
subroutine hist_in_0d(var, item, desc, unit)
Wrapper routine of HIST_reg+HIST_put 0D.
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
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]
real(rp), dimension(:), allocatable, public grid_ufz
face coordinate [m]: z, local=global
real(dp), public time_dtsec
time interval of model [sec]
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]
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
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
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
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
subroutine, public hist_setup
Setup.
subroutine hist_put_0d(itemid, var)
Put 1D data to history buffer.
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
logical, public io_lnml
output log or not? (for namelist, this process)
real(rp), dimension(:,:), allocatable, public real_lon
longitude [rad,0-2pi]
logical, public interp_available
vertical interpolation is available? (have meaning?)
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
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
integer, public io_fid_conf
Config file ID.
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.
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
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.
real(rp), dimension(:), allocatable, public grid_cbfxg
center buffer factor [0-1]: x, global
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