42 public :: file_cartesc_check_coordinates
44 public :: file_cartesc_get_size
50 public :: file_cartesc_write_var
51 public :: file_cartesc_read
52 public :: file_cartesc_write
56 interface file_cartesc_check_coordinates
57 module procedure file_cartesc_check_coordinates_name
58 module procedure file_cartesc_check_coordinates_id
59 end interface file_cartesc_check_coordinates
61 interface file_cartesc_get_size
63 module procedure file_cartesc_get_size_name
64 end interface file_cartesc_get_size
66 interface file_cartesc_read
68 module procedure file_cartesc_read_2d
69 module procedure file_cartesc_read_3d
70 module procedure file_cartesc_read_4d
71 module procedure file_cartesc_read_var_1d
77 end interface file_cartesc_read
79 interface file_cartesc_write
81 module procedure file_cartesc_write_2d
82 module procedure file_cartesc_write_3d
83 module procedure file_cartesc_write_3d_t
84 module procedure file_cartesc_write_4d
85 end interface file_cartesc_write
87 interface file_cartesc_write_var
93 end interface file_cartesc_write_var
100 integer :: size_global (1)
101 integer :: start_global(1)
102 integer :: halo_global (2)
103 integer :: halo_local (2)
119 real(
rp),
private :: file_cartesc_datacheck_criteria
122 character(len=H_SHORT) :: name
124 character(len=H_SHORT) :: dims(3)
127 character(len=H_SHORT) :: area
128 character(len=H_SHORT) :: area_x
129 character(len=H_SHORT) :: area_y
130 character(len=H_SHORT) :: volume
131 character(len=H_SHORT) :: location
132 character(len=H_SHORT) :: grid
134 integer,
parameter :: file_cartesc_ndims = 44
135 type(dims) :: file_cartesc_dims(file_cartesc_ndims)
140 real(
rp),
private,
allocatable :: axis_hgt (:,:,:)
141 real(
rp),
private,
allocatable :: axis_hgtwxy(:,:,:)
143 real(
rp),
private,
allocatable :: axis_lon (:,:)
144 real(
rp),
private,
allocatable :: axis_lonuy(:,:)
145 real(
rp),
private,
allocatable :: axis_lonxv(:,:)
146 real(
rp),
private,
allocatable :: axis_lonuv(:,:)
147 real(
rp),
private,
allocatable :: axis_lat (:,:)
148 real(
rp),
private,
allocatable :: axis_latuy(:,:)
149 real(
rp),
private,
allocatable :: axis_latxv(:,:)
150 real(
rp),
private,
allocatable :: axis_latuv(:,:)
152 real(
rp),
private,
allocatable :: axis_topo (:,:)
153 real(
rp),
private,
allocatable :: axis_lsmask(:,:)
155 real(
rp),
private,
allocatable :: axis_area (:,:)
156 real(
rp),
private,
allocatable :: axis_areazuy_x(:,:,:)
157 real(
rp),
private,
allocatable :: axis_areazxv_y(:,:,:)
158 real(
rp),
private,
allocatable :: axis_areawuy_x(:,:,:)
159 real(
rp),
private,
allocatable :: axis_areawxv_y(:,:,:)
160 real(
rp),
private,
allocatable :: axis_areauy (:,:)
161 real(
rp),
private,
allocatable :: axis_areazxy_x(:,:,:)
162 real(
rp),
private,
allocatable :: axis_areazuv_y(:,:,:)
163 real(
rp),
private,
allocatable :: axis_areaxv (:,:)
164 real(
rp),
private,
allocatable :: axis_areazuv_x(:,:,:)
165 real(
rp),
private,
allocatable :: axis_areazxy_y(:,:,:)
167 real(
rp),
private,
allocatable :: axis_vol (:,:,:)
168 real(
rp),
private,
allocatable :: axis_volwxy(:,:,:)
169 real(
rp),
private,
allocatable :: axis_volzuy(:,:,:)
170 real(
rp),
private,
allocatable :: axis_volzxv(:,:,:)
172 real(
rp),
private,
allocatable :: axis_volo(:,:,:)
173 real(
rp),
private,
allocatable :: axis_voll(:,:,:)
174 real(
rp),
private,
allocatable :: axis_volu(:,:,:)
181 integer,
private :: nfiles = 0
185 integer,
private,
target :: startxy (3), countxy (3)
186 integer,
private,
target :: startzx (2), countzx (2)
187 integer,
private,
target :: startzxy (4), countzxy (4)
188 integer,
private,
target :: startzhxy (4), countzhxy (4)
189 integer,
private,
target :: startocean(4), countocean(4)
190 integer,
private,
target :: startland (4), countland (4)
191 integer,
private,
target :: starturban(4), counturban(4)
193 integer,
private :: isb2, ieb2, jsb2, jeb2
196 integer,
private :: etype
199 integer,
private :: centertypexy
200 integer,
private :: centertypezx
201 integer,
private :: centertypezxy
202 integer,
private :: centertypezhxy
203 integer,
private :: centertypeocean
204 integer,
private :: centertypeland
205 integer,
private :: centertypeurban
207 logical,
private :: set_coordinates = .false.
209 logical,
private :: prof = .false.
226 namelist / param_file_cartesc / &
227 file_cartesc_datacheck_criteria
237 log_info(
"FILE_CARTESC_setup",*)
'Setup'
239 file_cartesc_datacheck_criteria = 0.1_rp**(
rp)
243 read(
io_fid_conf,nml=param_file_cartesc,iostat=ierr)
245 log_info(
"FILE_CARTESC_setup",*)
'Not found namelist. Default used.'
246 elseif( ierr > 0 )
then
247 log_error(
"FILE_CARTESC_setup",*)
'Not appropriate names in namelist PARAM_FILE_CARTESC. Check!'
250 log_nml(param_file_cartesc)
253 log_info(
"FILE_CARTESC_setup",*)
'NetCDF header information '
254 log_info_cont(*)
'Data source : ', trim(
h_source)
258 log_info(
"FILE_CARTESC_setup",*)
'Data consistency criteria : ', &
259 '(file-internal)/internal = ', file_cartesc_datacheck_criteria
274 allocate( axis_hgt(
kmax ,im,jm) )
275 allocate( axis_hgtwxy(
kmax+1,im,jm) )
277 allocate( axis_lon(im,jm) )
278 allocate( axis_lonuy(im,jm) )
279 allocate( axis_lonxv(im,jm) )
280 allocate( axis_lonuv(im,jm) )
281 allocate( axis_lat(im,jm) )
282 allocate( axis_latuy(im,jm) )
283 allocate( axis_latxv(im,jm) )
284 allocate( axis_latuv(im,jm) )
286 allocate( axis_topo(im,jm) )
287 allocate( axis_lsmask(im,jm) )
289 allocate( axis_area( im,jm) )
290 allocate( axis_areazuy_x(
kmax, im,jm) )
291 allocate( axis_areazxv_y(
kmax, im,jm) )
292 allocate( axis_areawuy_x(
kmax+1,im,jm) )
293 allocate( axis_areawxv_y(
kmax+1,im,jm) )
294 allocate( axis_areauy( im,jm) )
295 allocate( axis_areazxy_x(
kmax, im,jm) )
296 allocate( axis_areazuv_y(
kmax, im,jm) )
297 allocate( axis_areaxv( im,jm) )
298 allocate( axis_areazuv_x(
kmax, im,jm) )
299 allocate( axis_areazxy_y(
kmax, im,jm) )
301 allocate( axis_vol(
kmax ,im,jm) )
302 allocate( axis_volwxy(
kmax+1,im,jm) )
303 allocate( axis_volzuy(
kmax ,im,jm) )
304 allocate( axis_volzxv(
kmax ,im,jm) )
306 allocate( axis_volo(
okmax,im,jm) )
307 allocate( axis_voll(
lkmax,im,jm) )
308 allocate( axis_volu(
ukmax,im,jm) )
312 write_buf_amount(:) = 0
330 deallocate( axis_hgt )
331 deallocate( axis_hgtwxy )
333 deallocate( axis_lon )
334 deallocate( axis_lonuy )
335 deallocate( axis_lonxv )
336 deallocate( axis_lonuv )
337 deallocate( axis_lat )
338 deallocate( axis_latuy )
339 deallocate( axis_latxv )
340 deallocate( axis_latuv )
342 deallocate( axis_topo )
343 deallocate( axis_lsmask )
345 deallocate( axis_area )
346 deallocate( axis_areazuy_x )
347 deallocate( axis_areazxv_y )
348 deallocate( axis_areawuy_x )
349 deallocate( axis_areawxv_y )
350 deallocate( axis_areauy )
351 deallocate( axis_areazxy_x )
352 deallocate( axis_areazuv_y )
353 deallocate( axis_areaxv )
354 deallocate( axis_areazuv_x )
355 deallocate( axis_areazxy_y )
357 deallocate( axis_vol )
358 deallocate( axis_volwxy )
359 deallocate( axis_volzuy )
360 deallocate( axis_volzxv )
362 deallocate( axis_volo )
363 deallocate( axis_voll )
364 deallocate( axis_volu )
368 set_coordinates = .false.
378 subroutine file_cartesc_get_size_name( &
380 KMAX, OKMAX, LKMAX, UKMAX, &
382 KHALO, IHALO, JHALO, &
386 character(len=*),
intent(in) :: basename
392 logical,
intent(in),
optional :: aggregate
398 aggregate=aggregate )
406 end subroutine file_cartesc_get_size_name
409 KMAX, OKMAX, LKMAX, UKMAX, &
411 KHALO, IHALO, JHALO )
415 integer,
intent(in) :: fid
417 integer,
intent(out) :: KMAX, OKMAX, LKMAX, UKMAX
418 integer,
intent(out) :: IMAXG, JMAXG
419 integer,
intent(out) :: KHALO, IHALO, JHALO
424 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_kmax", buf(:) )
426 call file_get_attribute( fid,
"global",
"scale_ocean_grid_cartesC_index_kmax", buf(:), existed=existed )
432 call file_get_attribute( fid,
"global",
"scale_land_grid_cartesC_index_kmax", buf(:), existed=existed )
438 call file_get_attribute( fid,
"global",
"scale_urban_grid_cartesC_index_kmax", buf(:), existed=existed )
445 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_imaxg", buf(:) )
447 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_jmaxg", buf(:) )
450 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_khalo", buf(:) )
452 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_ihalo", buf(:) )
454 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_jhalo", buf(:) )
464 LON, LONUY, LONXV, LONUV, &
465 LAT, LATUY, LATXV, LATUV, &
471 real(
rp),
intent(in) :: cz(
ka,
ia,
ja)
472 real(
rp),
intent(in) :: fz(0:
ka,
ia,
ja)
473 real(
rp),
intent(in) :: lon (
ia,
ja)
474 real(
rp),
intent(in) :: lonuy(0:
ia,
ja)
475 real(
rp),
intent(in) :: lonxv(
ia,0:
ja)
476 real(
rp),
intent(in) :: lonuv(0:
ia,0:
ja)
477 real(
rp),
intent(in) :: lat (
ia,
ja)
478 real(
rp),
intent(in) :: latuy(0:
ia,
ja)
479 real(
rp),
intent(in) :: latxv(
ia,0:
ja)
480 real(
rp),
intent(in) :: latuv(0:
ia,0:
ja)
481 real(
rp),
intent(in) :: topo (
ia,
ja)
482 real(
rp),
intent(in) :: lsmask(
ia,
ja)
485 axis_hgt(:,:,:) = cz(
ks :
ke,isb2:ieb2,jsb2:jeb2)
486 axis_hgtwxy(:,:,:) = fz(
ks-1:
ke,isb2:ieb2,jsb2:jeb2)
488 axis_lon(:,:) = lon(isb2:ieb2,jsb2:jeb2) / d2r
489 axis_lonuy(:,:) = lonuy(isb2:ieb2,jsb2:jeb2) / d2r
490 axis_lonxv(:,:) = lonxv(isb2:ieb2,jsb2:jeb2) / d2r
491 axis_lonuv(:,:) = lonuv(isb2:ieb2,jsb2:jeb2) / d2r
492 axis_lat(:,:) = lat(isb2:ieb2,jsb2:jeb2) / d2r
493 axis_latuy(:,:) = latuy(isb2:ieb2,jsb2:jeb2) / d2r
494 axis_latxv(:,:) = latxv(isb2:ieb2,jsb2:jeb2) / d2r
495 axis_latuv(:,:) = latuv(isb2:ieb2,jsb2:jeb2) / d2r
497 axis_topo(:,:) = topo(isb2:ieb2,jsb2:jeb2)
498 axis_lsmask(:,:) = lsmask(isb2:ieb2,jsb2:jeb2)
500 set_coordinates = .true.
508 AREA, AREAZUY_X, AREAZXV_Y, &
509 AREAWUY_X, AREAWXV_Y, &
510 AREAUY, AREAZXY_X, AREAZUV_Y, &
511 AREAXV, AREAZUV_X, AREAZXY_Y, &
512 VOL, VOLWXY, VOLZUY, VOLZXV )
516 real(
rp),
intent(in) :: area (
ia,
ja)
517 real(
rp),
intent(in) :: areazuy_x(
ka,
ia,
ja)
518 real(
rp),
intent(in) :: areazxv_y(
ka,
ia,
ja)
519 real(
rp),
intent(in) :: areawuy_x(0:
ka,
ia,
ja)
520 real(
rp),
intent(in) :: areawxv_y(0:
ka,
ia,
ja)
521 real(
rp),
intent(in) :: areauy (
ia,
ja)
522 real(
rp),
intent(in) :: areazxy_x(
ka,
ia,
ja)
523 real(
rp),
intent(in) :: areazuv_y(
ka,
ia,
ja)
524 real(
rp),
intent(in) :: areaxv (
ia,
ja)
525 real(
rp),
intent(in) :: areazuv_x(
ka,
ia,
ja)
526 real(
rp),
intent(in) :: areazxy_y(
ka,
ia,
ja)
527 real(
rp),
intent(in) :: vol (
ka,
ia,
ja)
528 real(
rp),
intent(in) :: volwxy(0:
ka,
ia,
ja)
529 real(
rp),
intent(in) :: volzuy(
ka,
ia,
ja)
530 real(
rp),
intent(in) :: volzxv(
ka,
ia,
ja)
532 axis_area(:,:) = area( isb2:ieb2,jsb2:jeb2)
533 axis_areazuy_x(:,:,:) = areazuy_x(
ks :
ke,isb2:ieb2,jsb2:jeb2)
534 axis_areazxv_y(:,:,:) = areazxv_y(
ks :
ke,isb2:ieb2,jsb2:jeb2)
535 axis_areawuy_x(:,:,:) = areawuy_x(
ks-1:
ke,isb2:ieb2,jsb2:jeb2)
536 axis_areawxv_y(:,:,:) = areawxv_y(
ks-1:
ke,isb2:ieb2,jsb2:jeb2)
537 axis_areauy(:,:) = areauy( isb2:ieb2,jsb2:jeb2)
538 axis_areazxy_x(:,:,:) = areazxy_x(
ks :
ke,isb2:ieb2,jsb2:jeb2)
539 axis_areazuv_y(:,:,:) = areazuv_y(
ks :
ke,isb2:ieb2,jsb2:jeb2)
540 axis_areaxv(:,:) = areaxv( isb2:ieb2,jsb2:jeb2)
541 axis_areazuv_x(:,:,:) = areazuv_x(
ks :
ke,isb2:ieb2,jsb2:jeb2)
542 axis_areazxy_y(:,:,:) = areazxy_y(
ks :
ke,isb2:ieb2,jsb2:jeb2)
544 axis_vol(:,:,:) = vol(
ks :
ke,isb2:ieb2,jsb2:jeb2)
545 axis_volwxy(:,:,:) = volwxy(
ks-1:
ke,isb2:ieb2,jsb2:jeb2)
546 axis_volzuy(:,:,:) = volzuy(
ks :
ke,isb2:ieb2,jsb2:jeb2)
547 axis_volzxv(:,:,:) = volzxv(
ks :
ke,isb2:ieb2,jsb2:jeb2)
561 axis_volo(:,:,:) = vol(
oks:
oke,isb2:ieb2,jsb2:jeb2)
575 axis_voll(:,:,:) = vol(
lks:
lke,isb2:ieb2,jsb2:jeb2)
589 axis_volu(:,:,:) = vol(
uks:
uke,isb2:ieb2,jsb2:jeb2)
596 subroutine file_cartesc_check_coordinates_name( &
598 atmos, ocean, land, urban, &
602 character(len=*),
intent(in) :: basename
603 logical,
intent(in),
optional :: atmos
604 logical,
intent(in),
optional :: ocean
605 logical,
intent(in),
optional :: land
606 logical,
intent(in),
optional :: urban
607 logical,
intent(in),
optional :: transpose
613 logical :: transpose_
624 if(
present(atmos) ) atmos_ = atmos
625 if(
present(ocean) ) ocean_ = ocean
626 if(
present(land ) ) land_ = land
627 if(
present(urban) ) urban_ = urban
628 if(
present(transpose) ) transpose_ = transpose
633 call file_cartesc_check_coordinates_id( fid, &
634 atmos_, ocean_, land_, urban_, &
638 end subroutine file_cartesc_check_coordinates_name
642 subroutine file_cartesc_check_coordinates_id( &
644 atmos, ocean, land, urban, &
658 integer,
intent(in) :: fid
659 logical,
intent(in),
optional :: atmos
660 logical,
intent(in),
optional :: ocean
661 logical,
intent(in),
optional :: land
662 logical,
intent(in),
optional :: urban
663 logical,
intent(in),
optional :: transpose
669 logical :: transpose_
671 real(
rp) :: buffer_z (
ka)
672 real(
rp) :: buffer_x (
ia)
673 real(
rp) :: buffer_y (
ja)
674 real(
rp) :: buffer_xy (
ia,
ja)
680 integer :: xsb, xeb, ysb, yeb
684 log_info(
"FILE_CARTESC_check_coordinates_id",*)
'Check consistency of axis '
692 if(
present(atmos) ) atmos_ = atmos
693 if(
present(ocean) ) ocean_ = ocean
694 if(
present(land ) ) land_ = land
695 if(
present(urban) ) urban_ = urban
696 if(
present(transpose) ) transpose_ = transpose
704 call file_cartesc_read_var_1d( fid,
'x',
'X', buffer_x(:) )
705 call file_cartesc_read_var_1d( fid,
'y',
'Y', buffer_y(:) )
710 if ( set_coordinates )
then
713 call check_2d( axis_lon(xsb:xeb,ysb:yeb), buffer_xy(
isb:
ieb,
jsb:
jeb),
'lon' )
717 call check_2d( axis_lat(xsb:xeb,ysb:yeb), buffer_xy(
isb:
ieb,
jsb:
jeb),
'lat' )
721 call file_cartesc_read_var_1d( fid,
'z',
'Z', buffer_z(:) )
722 if ( .not. transpose_ )
then
727 if ( .not. transpose_ )
then
728 call check_3d( axis_hgt(:,xsb:xeb,ysb:yeb), buffer_zxy(
ks:
ke,
isb:
ieb,
jsb:
jeb),
'height', transpose_ )
733 call file_cartesc_read_var_1d( fid,
'oz',
'OZ', buffer_o(:) )
739 call file_cartesc_read_var_1d( fid,
'lz',
'LZ', buffer_l(:) )
745 call file_cartesc_read_var_1d( fid,
'uz',
'UZ', buffer_u(:) )
751 end subroutine file_cartesc_check_coordinates_id
769 character(len=*),
intent(in) :: basename
770 integer,
intent(out) :: fid
771 logical,
intent(in),
optional :: single
772 logical,
intent(in),
optional :: aggregate
780 aggregate=aggregate, &
791 basename, title, datatype, &
795 append, aggregate, single )
819 character(len=*),
intent(in) :: basename
820 character(len=*),
intent(in) :: title
821 character(len=*),
intent(in) :: datatype
822 integer,
intent(out) :: fid
823 integer,
intent(in),
optional :: date(6)
824 real(
dp),
intent(in),
optional :: subsec
825 logical,
intent(in),
optional :: append
826 logical,
intent(in),
optional :: haszcoord
827 logical,
intent(in),
optional :: aggregate
828 logical,
intent(in),
optional :: single
832 character(len=34) :: tunits
833 character(len=H_SHORT) :: calendar
835 integer :: rank_x, rank_y
836 integer :: num_x, num_y
837 logical :: fileexisted
838 logical :: aggregate_
846 if (
present(single) )
then
853 if ( datatype ==
'REAL8' )
then
855 elseif( datatype ==
'REAL4' )
then
860 elseif(
rp == 4 )
then
863 log_error(
"FILE_CARTESC_create",*)
'unsupported data type. Check!', trim(datatype)
869 if (
present(append) )
then
874 if (
present(date) )
then
877 date_(:) = nowdate(:)
879 if ( date_(1) > 0 )
then
889 if (
present(aggregate) )
then
890 aggregate_ = aggregate
903 aggregate = aggregate_, &
904 time_units = tunits, &
905 calendar = calendar, &
909 if ( fid > 0 .and. (.not. fileexisted) )
then
914 file_axes_written(fid) = .false.
916 if (
present( haszcoord ) )
then
917 file_haszcoord(fid) = haszcoord
919 file_haszcoord(fid) = .true.
922 if ( aggregate_ )
then
934 if (
present( subsec ) )
then
947 subsec_, tunits, calendar )
951 file_haszcoord(fid) )
973 integer,
intent(in) :: fid
985 if ( .NOT. file_axes_written(fid) )
then
998 file_haszcoord(fid), &
1001 file_axes_written(fid) = .true.
1025 integer,
intent(in) :: fid
1053 integer,
intent(in) :: fid
1062 if ( write_buf_amount(fid) > 0 )
then
1064 write_buf_amount(fid) = 0
1082 basename, varname, &
1089 character(len=*),
intent(in) :: basename
1090 character(len=*),
intent(in) :: varname
1091 character(len=*),
intent(in) :: dim_type
1093 real(RP),
intent(out) :: var(:)
1095 integer,
intent(in),
optional :: step
1096 logical,
intent(in),
optional :: aggregate
1097 logical,
intent(in),
optional :: allow_missing
1104 aggregate=aggregate )
1106 call file_cartesc_read_var_1d( fid, varname, dim_type, &
1109 allow_missing=allow_missing )
1118 subroutine file_cartesc_read_2d( &
1119 basename, varname, &
1126 character(len=*),
intent(in) :: basename
1127 character(len=*),
intent(in) :: varname
1128 character(len=*),
intent(in) :: dim_type
1130 real(RP),
intent(out) :: var(:,:)
1132 integer,
intent(in),
optional :: step
1133 logical,
intent(in),
optional :: aggregate
1134 logical,
intent(in),
optional :: allow_missing
1141 aggregate=aggregate )
1146 allow_missing=allow_missing )
1151 end subroutine file_cartesc_read_2d
1155 subroutine file_cartesc_read_3d( &
1156 basename, varname, &
1163 character(len=*),
intent(in) :: basename
1164 character(len=*),
intent(in) :: varname
1165 character(len=*),
intent(in) :: dim_type
1167 real(RP),
intent(out) :: var(:,:,:)
1169 integer,
intent(in),
optional :: step
1170 logical,
intent(in),
optional :: aggregate
1171 logical,
intent(in),
optional :: allow_missing
1178 aggregate=aggregate )
1183 allow_missing=allow_missing )
1188 end subroutine file_cartesc_read_3d
1192 subroutine file_cartesc_read_4d( &
1193 basename, varname, &
1199 character(len=*),
intent(in) :: basename
1200 character(len=*),
intent(in) :: varname
1201 character(len=*),
intent(in) :: dim_type
1202 integer,
intent(in) :: step
1204 real(RP),
intent(out) :: var(:,:,:,:)
1206 logical,
intent(in),
optional :: aggregate
1207 logical,
intent(in),
optional :: allow_missing
1214 aggregate=aggregate )
1218 allow_missing=allow_missing )
1223 end subroutine file_cartesc_read_4d
1227 subroutine file_cartesc_read_var_1d( &
1245 integer,
intent(in) :: fid
1246 character(len=*),
intent(in) :: varname
1247 character(len=*),
intent(in) :: dim_type
1249 real(RP),
intent(out) :: var(:)
1251 integer,
intent(in),
optional :: step
1252 logical,
intent(in),
optional :: allow_missing
1255 integer :: dim1_S, dim1_E
1264 log_info(
"FILE_CARTESC_read_var_1D",
'(1x,2A)')
'Read from file (1D), name : ', trim(varname)
1268 if ( dim_type ==
'Z' )
then
1273 elseif( dim_type ==
'OZ' )
then
1278 elseif( dim_type ==
'LZ' )
then
1283 elseif( dim_type ==
'UZ' )
then
1288 elseif( dim_type ==
'X' .OR. dim_type ==
'CX' )
then
1293 elseif( dim_type ==
'Y' .OR. dim_type ==
'CY' )
then
1299 log_error(
"FILE_CARTESC_read_var_1D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1303 if (
size(var) .ne. vsize )
then
1304 log_error(
"FILE_CARTESC_read_var_1D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1307 count(1) = dim1_e - dim1_s + 1
1308 call file_read( fid, varname, &
1309 var(dim1_s:dim1_e), &
1310 step=step, allow_missing=allow_missing, &
1311 ntypes=count(1), dtype=etype, start=start, count=count )
1314 if ( dim_type ==
'Z' )
then
1318 elseif( dim_type ==
'OZ' )
then
1322 elseif( dim_type ==
'LZ' )
then
1326 elseif( dim_type ==
'UZ' )
then
1330 elseif( dim_type ==
'X' )
then
1334 elseif( dim_type ==
'CX' )
then
1338 elseif( dim_type ==
'Y' )
then
1342 elseif( dim_type ==
'CY' )
then
1347 log_error(
"FILE_CARTESC_read_var_1D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1351 if (
size(var) .ne. vsize )
then
1352 log_error(
"FILE_CARTESC_read_var_1D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1355 call file_read( fid, varname, var(dim1_s:dim1_e), step=step )
1361 end subroutine file_cartesc_read_var_1d
1379 integer,
intent(in) :: fid
1380 character(len=*),
intent(in) :: varname
1381 character(len=*),
intent(in) :: dim_type
1383 real(RP),
intent(out) :: var(:,:)
1385 integer,
intent(in),
optional :: step
1386 logical,
intent(in),
optional :: allow_missing
1389 integer :: ntypes, dtype
1390 integer,
pointer :: start(:), count(:)
1391 integer :: dim1_S, dim1_E
1392 integer :: dim2_S, dim2_E
1394 real(RP),
allocatable :: buf(:,:)
1401 log_info(
"FILE_CARTESC_read_var_2D",
'(1x,2A)')
'Read from file (2D), name : ', trim(varname)
1406 if ( dim_type ==
'XY' )
then
1412 elseif( dim_type ==
'ZX' )
then
1417 dtype = centertypezx
1421 log_error(
"FILE_CARTESC_read_var_2D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1425 if (
size(var) .ne. vsize )
then
1426 log_error(
"FILE_CARTESC_read_var_2D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1429 call file_read( fid, varname, &
1431 step=step, allow_missing=allow_missing, &
1432 ntypes=ntypes, dtype=dtype, start=start, count=count )
1436 if ( dim_type ==
'XY' )
then
1442 elseif( dim_type ==
'ZX' )
then
1449 log_error(
"FILE_CARTESC_read_var_2D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1453 if (
size(var) .ne. vsize )
then
1454 log_error(
"FILE_CARTESC_read_var_2D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1457 allocate( buf(dim1_s:dim1_e,dim2_s:dim2_e) )
1458 call file_read( fid, varname, buf(:,:), step=step )
1460 var(dim1_s:dim1_e,dim2_s:dim2_e) = buf(:,:)
1491 integer,
intent(in) :: fid
1492 character(len=*),
intent(in) :: varname
1493 character(len=*),
intent(in) :: dim_type
1495 real(RP),
intent(out) :: var(:,:,:)
1497 integer,
intent(in),
optional :: step
1498 logical,
intent(in),
optional :: allow_missing
1501 integer :: ntypes, dtype
1502 integer,
pointer :: start(:), count(:)
1503 integer :: dim1_S, dim1_E
1504 integer :: dim2_S, dim2_E
1505 integer :: dim3_S, dim3_E
1506 real(RP),
allocatable :: buf(:,:,:)
1513 log_info(
"FILE_CARTESC_read_var_3D",
'(1x,2A)')
'Read from file (3D), name : ', trim(varname)
1520 if( dim_type ==
'ZXY' &
1521 .or. dim_type ==
'ZXHY' &
1522 .or. dim_type ==
'ZXYH' )
then
1525 dtype = centertypezxy
1528 elseif( dim_type ==
'ZHXY' )
then
1531 dtype = centertypezhxy
1534 elseif( dim_type ==
'XYT' )
then
1535 if ( .not.
present(step) )
then
1536 log_error(
"FILE_CARTESC_read_var_3D",*)
'step is necessary for "XYT"'
1539 vsize =
ia *
ja * step
1540 ntypes =
ia *
ja * step
1546 elseif( dim_type ==
'OXY' )
then
1549 dtype = centertypeocean
1552 elseif( dim_type ==
'LXY' )
then
1555 dtype = centertypeland
1558 elseif( dim_type ==
'UXY' )
then
1561 dtype = centertypeurban
1565 log_error(
"FILE_CARTESC_read_var_3D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1569 if (
size(var) .ne. vsize )
then
1570 log_error(
"FILE_CARTESC_read_var_3D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1573 call file_read( fid, varname, &
1575 step=step, allow_missing=allow_missing, &
1576 ntypes=ntypes, dtype=dtype, start=start, count=count )
1579 if( dim_type ==
'ZXY' &
1580 .or. dim_type ==
'ZXHY' &
1581 .or. dim_type ==
'ZXYH' )
then
1589 elseif( dim_type ==
'ZHXY' )
then
1597 elseif( dim_type ==
'XYT' )
then
1598 if ( .not.
present(step) )
then
1599 log_error(
"FILE_CARTESC_read_var_3D",*)
'step is necessary for "XYT"'
1602 vsize =
ia *
ja * step
1609 elseif( dim_type ==
'OXY' )
then
1617 elseif( dim_type ==
'LXY' )
then
1625 elseif( dim_type ==
'UXY' )
then
1634 log_error(
"FILE_CARTESC_read_var_3D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1638 if (
size(var) .ne. vsize )
then
1639 log_error(
"FILE_CARTESC_read_var_3D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1643 allocate( buf(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e) )
1644 call file_read( fid, varname, buf(:,:,:), &
1645 step=step, allow_missing=allow_missing )
1647 var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e) = buf(:,:,:)
1678 integer,
intent(in) :: fid
1679 character(len=*),
intent(in) :: varname
1680 character(len=*),
intent(in) :: dim_type
1681 integer,
intent(in) :: step
1683 real(RP),
intent(out) :: var(:,:,:,:)
1685 logical,
intent(in),
optional :: allow_missing
1689 integer,
pointer :: start(:), count(:)
1690 integer :: dim1_S, dim1_E
1691 integer :: dim2_S, dim2_E
1692 integer :: dim3_S, dim3_E
1693 integer :: dim4_S, dim4_E
1695 real(RP),
allocatable :: buf(:,:,:,:)
1702 log_info(
"FILE_CARTESC_read_var_4D",
'(1x,2A)')
'Read from file (4D), name : ', trim(varname)
1706 if ( dim_type ==
'ZXYT' &
1707 .or. dim_type ==
'ZXHYT' &
1708 .or. dim_type ==
'ZXYHT' )
then
1709 vsize =
ka *
ia *
ja * step
1710 dtype = centertypezxy
1713 elseif ( dim_type ==
'ZHXYT' )
then
1714 vsize =
ka *
ia *
ja * step
1715 dtype = centertypezhxy
1718 elseif ( dim_type ==
'OXYT' )
then
1720 dtype = centertypeocean
1723 elseif ( dim_type ==
'LXYT' )
then
1725 dtype = centertypeland
1728 elseif ( dim_type ==
'UXYT' )
then
1730 dtype = centertypeurban
1734 log_error(
"FILE_CARTESC_read_var_4D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1738 if (
size(var) .ne. vsize )
then
1739 log_error(
"FILE_CARTESC_read_var_4D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1744 call file_read( fid, varname, &
1746 allow_missing=allow_missing, &
1747 ntypes=step, dtype=dtype, start=start, count=count )
1750 if ( dim_type ==
'ZXYT' &
1751 .or. dim_type ==
'ZXHYT' &
1752 .or. dim_type ==
'ZXYHT' )
then
1753 vsize =
ka *
ia *
ja * step
1760 elseif ( dim_type ==
'ZHXYT' )
then
1761 vsize =
ka *
ia *
ja * step
1768 elseif ( dim_type ==
'OXYT' )
then
1776 elseif ( dim_type ==
'LXYT' )
then
1784 elseif ( dim_type ==
'UXYT' )
then
1793 log_error(
"FILE_CARTESC_read_var_4D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
1797 if (
size(var) .ne. vsize )
then
1798 log_error(
"FILE_CARTESC_read_var_4D",*)
'size of var is invalid: ', trim(varname),
size(var), vsize
1804 allocate( buf(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e,dim4_s:dim4_e) )
1805 call file_read( fid, varname, &
1808 var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e,dim4_s:dim4_e) = buf(:,:,:,:)
1830 file_get_datainfo, &
1831 file_get_attribute, &
1836 integer,
intent(in) :: fid
1837 character(len=*),
intent(in) :: varname
1839 real(RP),
intent(out) :: var(:,:)
1841 integer,
intent(in),
optional :: step
1843 logical,
intent(out),
optional :: existed
1849 character(len=H_SHORT) :: dnames(2)
1863 log_info(
"FILE_CARTESC_read_auto_2D",
'(1x,2A)')
'Read from file (2D), name : ', trim(varname)
1865 call file_get_datainfo( fid, varname, dim_name=dnames(:), existed=existed2 )
1867 if (
present( existed ) )
then
1869 if ( .not. existed2 )
then
1875 if ( .not. existed2 )
then
1876 log_error(
"FILE_CARTESC_read_auto_2D",*)
'variable not found: ', trim(varname)
1880 call file_get_shape( fid, varname, dims(:) )
1884 if ( nx==dims(1) .and. ny==dims(2) )
then
1888 call file_get_attribute( fid, dnames(n),
"halo_local", halos(:), existed=existed2 )
1889 if ( existed2 )
then
1890 start(n) = halos(1) + 1
1896 count(:) = (/nx,ny/)
1898 call file_read( fid, varname, var(:,:), step=step, start=start(:), count=count(:) )
1917 file_get_datainfo, &
1918 file_get_attribute, &
1923 integer,
intent(in) :: fid
1924 character(len=*),
intent(in) :: varname
1926 real(RP),
intent(out) :: var(:,:,:)
1928 integer,
intent(in),
optional :: step
1930 logical,
intent(out),
optional :: existed
1936 character(len=H_SHORT) :: dnames(3)
1940 real(RP),
allocatable :: buf(:,:,:)
1941 integer :: nx, ny, nz
1952 log_info(
"FILE_CARTESC_read_auto_3D",
'(1x,2A)')
'Read from file (3D), name : ', trim(varname)
1954 call file_get_datainfo( fid, varname, dim_name=dnames(:), existed=existed2 )
1956 if (
present(existed) )
then
1958 if ( .not. existed2 )
then
1964 if ( .not. existed2 )
then
1965 log_error(
"FILE_CARTESC_read_auto_3D",*)
'variable not found: ', trim(varname)
1969 call file_get_shape( fid, varname, dims(:) )
1974 if ( ( dnames(1)(1:1)==
"z" .or. dnames(1)(2:2)==
"z" ) .and. dnames(2)(1:1)==
"x" .and. dnames(3)(1:1)==
"y" )
then
1975 if ( nz==dims(1) .and. nx==dims(2) .and. ny==dims(3) )
then
1976 start(:) = (/1,1,1/)
1977 else if ( dnames(1)==
"zh" .and. nz+1==dims(1) .and. nx==dims(2) .and. ny==dims(3) )
then
1978 start(:) = (/2,1,1/)
1981 call file_get_attribute( fid, dnames(n),
"halo_local", halos(:), existed=existed2 )
1982 if ( existed2 )
then
1983 start(n) = halos(1) + 1
1984 else if ( dnames(n)==
"zh" )
then
1991 count(:) = (/nz,nx,ny/)
1992 call file_read( fid, varname, var(:,:,:), step=step, start=start(:), count=count(:) )
1995 else if ( dnames(1)(1:1)==
"x" .and. dnames(2)(1:1)==
"y" .and. ( dnames(3)(1:1)==
"z" .or. dnames(3)(2:2)==
"z" ) )
then
1996 allocate( buf(nx,ny,nz) )
1997 if ( nx==dims(1) .and. ny==dims(2) .and. nz==dims(3) )
then
1998 start(:) = (/1,1,1/)
1999 else if ( nx==dims(1) .and. ny==dims(2) .and. nz+1==dims(3) .and. dnames(3)==
"zh" )
then
2000 start(:) = (/1,1,2/)
2003 call file_get_attribute( fid, dnames(n),
"halo_local", halos(:), existed=existed2 )
2004 if ( existed2 )
then
2005 start(n) = halos(1) + 1
2006 else if ( dnames(n)==
"zh" )
then
2013 count(:) = (/nx,ny,nz/)
2014 call file_read( fid, varname, buf(:,:,:), step=step, start=start(:), count=count(:) )
2022 var(k,i,j) = buf(i,j,k)
2029 log_error(
"FILE_CARTESC_read_auto_3D",*)
'invalid dimension'
2047 varname, desc, unit, &
2048 dim_type, datatype, &
2050 append, aggregate, &
2055 real(RP),
intent(in) :: var(:)
2056 character(len=*),
intent(in) :: basename
2057 character(len=*),
intent(in) :: title
2058 character(len=*),
intent(in) :: varname
2059 character(len=*),
intent(in) :: desc
2060 character(len=*),
intent(in) :: unit
2061 character(len=*),
intent(in) :: dim_type
2062 character(len=*),
intent(in) :: datatype
2064 integer,
intent(in),
optional :: date(6)
2065 real(DP),
intent(in),
optional :: subsec
2066 logical,
intent(in),
optional :: append
2067 logical,
intent(in),
optional :: aggregate
2068 character(len=*),
intent(in),
optional :: standard_name
2069 character(len=*),
intent(in),
optional :: cell_measures
2074 log_info(
"FILE_CARTESC_write_1D",
'(1x,2A)')
'Write to file (1D), name : ', trim(varname)
2078 date=date, subsec=subsec, &
2079 append=append, aggregate=aggregate, single=.true. )
2083 standard_name=standard_name, &
2084 cell_measures=cell_measures )
2095 subroutine file_cartesc_write_2d( &
2098 varname, desc, unit, &
2099 dim_type, datatype, &
2101 fill_halo, haszcoord, &
2102 append, aggregate, &
2107 real(RP),
intent(in) :: var(:,:)
2108 character(len=*),
intent(in) :: basename
2109 character(len=*),
intent(in) :: title
2110 character(len=*),
intent(in) :: varname
2111 character(len=*),
intent(in) :: desc
2112 character(len=*),
intent(in) :: unit
2113 character(len=*),
intent(in) :: dim_type
2114 character(len=*),
intent(in) :: datatype
2116 integer,
intent(in),
optional :: date(6)
2117 real(DP),
intent(in),
optional :: subsec
2118 logical,
intent(in),
optional :: fill_halo
2119 logical,
intent(in),
optional :: haszcoord
2120 logical,
intent(in),
optional :: append
2121 logical,
intent(in),
optional :: aggregate
2122 character(len=*),
intent(in),
optional :: standard_name
2123 character(len=*),
intent(in),
optional :: cell_measures
2128 log_info(
"FILE_CARTESC_write_2D",
'(1x,2A)')
'Write to file (2D), name : ', trim(varname)
2132 date=date, subsec=subsec, &
2133 haszcoord=haszcoord, &
2134 append=append, aggregate=aggregate )
2138 standard_name=standard_name, &
2139 cell_measures=cell_measures )
2146 end subroutine file_cartesc_write_2d
2150 subroutine file_cartesc_write_3d( &
2153 varname, desc, unit, &
2154 dim_type, datatype, &
2157 append, aggregate, &
2162 real(RP),
intent(in) :: var(:,:,:)
2163 character(len=*),
intent(in) :: basename
2164 character(len=*),
intent(in) :: title
2165 character(len=*),
intent(in) :: varname
2166 character(len=*),
intent(in) :: desc
2167 character(len=*),
intent(in) :: unit
2168 character(len=*),
intent(in) :: dim_type
2169 character(len=*),
intent(in) :: datatype
2171 integer,
intent(in),
optional :: date(6)
2172 real(DP),
intent(in),
optional :: subsec
2173 logical,
intent(in),
optional :: fill_halo
2174 logical,
intent(in),
optional :: append
2175 logical,
intent(in),
optional :: aggregate
2176 character(len=*),
intent(in),
optional :: standard_name
2177 character(len=*),
intent(in),
optional :: cell_measures
2182 log_info(
"FILE_CARTESC_write_3D",
'(1x,2A)')
'Write to file (3D), name : ', trim(varname)
2186 date=date, subsec=subsec, &
2187 append=append, aggregate=aggregate )
2191 standard_name=standard_name, &
2192 cell_measures=cell_measures )
2200 end subroutine file_cartesc_write_3d
2204 subroutine file_cartesc_write_3d_t( &
2207 varname, desc, unit, &
2208 dim_type, datatype, &
2210 timetarg, timeofs, &
2212 append, aggregate, &
2217 real(RP),
intent(in) :: var(:,:,:)
2218 character(len=*),
intent(in) :: basename
2219 character(len=*),
intent(in) :: title
2220 character(len=*),
intent(in) :: varname
2221 character(len=*),
intent(in) :: desc
2222 character(len=*),
intent(in) :: unit
2223 character(len=*),
intent(in) :: dim_type
2224 character(len=*),
intent(in) :: datatype
2225 real(DP),
intent(in) :: timeintv
2226 integer ,
intent(in) :: tsince(6)
2228 integer,
intent(in),
optional :: timetarg
2229 real(DP),
intent(in),
optional :: timeofs
2230 logical,
intent(in),
optional :: fill_halo
2231 logical,
intent(in),
optional :: append
2232 logical,
intent(in),
optional :: aggregate
2233 character(len=*),
intent(in),
optional :: standard_name
2234 character(len=*),
intent(in),
optional :: cell_measures
2242 log_info(
"FILE_CARTESC_write_3D_t",
'(1x,3A)')
'Write to file (3D), name : ', trim(varname),
'with time dimension'
2247 append=append, aggregate=aggregate )
2249 if (
present(timetarg) )
then
2252 nsteps =
size(var,3)
2256 standard_name=standard_name, &
2257 cell_measures=cell_measures, &
2258 timeintv=timeintv, nsteps=nsteps )
2263 timetarg, timeofs, fill_halo )
2266 end subroutine file_cartesc_write_3d_t
2270 subroutine file_cartesc_write_4d( &
2273 varname, desc, unit, &
2274 dim_type, datatype, &
2276 timetarg, timeofs, &
2278 append, aggregate, &
2283 real(RP),
intent(in) :: var(:,:,:,:)
2284 character(len=*),
intent(in) :: basename
2285 character(len=*),
intent(in) :: title
2286 character(len=*),
intent(in) :: varname
2287 character(len=*),
intent(in) :: desc
2288 character(len=*),
intent(in) :: unit
2289 character(len=*),
intent(in) :: dim_type
2290 character(len=*),
intent(in) :: datatype
2291 real(DP),
intent(in) :: timeintv
2292 integer,
intent(in) :: tsince(6)
2294 integer,
intent(in),
optional :: timetarg
2295 real(DP),
intent(in),
optional :: timeofs
2296 logical,
intent(in),
optional :: fill_halo
2297 logical,
intent(in),
optional :: append
2298 logical,
intent(in),
optional :: aggregate
2299 character(len=*),
intent(in),
optional :: standard_name
2300 character(len=*),
intent(in),
optional :: cell_measures
2308 log_info(
"FILE_CARTESC_write_4D",
'(1x,2A)')
'Write to file (4D), name : ', trim(varname)
2313 append=append, aggregate=aggregate )
2315 if (
present(timetarg) )
then
2318 nsteps =
size(var,3)
2322 standard_name=standard_name, &
2323 cell_measures=cell_measures, &
2324 timeintv=timeintv, nsteps=nsteps )
2329 timetarg, timeofs, fill_halo )
2332 end subroutine file_cartesc_write_4d
2338 prc_rank_x, prc_rank_y, &
2339 prc_num_x, prc_num_y, &
2340 prc_periodic_x, prc_periodic_y, &
2341 kmax, okmax, lkmax, ukmax, &
2343 khalo, ihalo, jhalo, &
2344 time, tunits, calendar )
2352 integer,
intent(in) :: fid
2353 integer,
intent(in) :: prc_rank_x, prc_rank_y
2354 integer,
intent(in) :: prc_num_x, prc_num_y
2355 logical,
intent(in) :: prc_periodic_x, prc_periodic_y
2359 real(
dp),
intent(in) :: time
2360 character(len=*),
intent(in) :: tunits
2361 character(len=*),
intent(in) :: calendar
2368 call file_set_attribute( fid,
"global",
"Conventions",
"CF-1.6" )
2372 call file_set_attribute( fid,
"global",
"scale_cartesC_prc_rank_x", (/prc_rank_x/) )
2373 call file_set_attribute( fid,
"global",
"scale_cartesC_prc_rank_y", (/prc_rank_y/) )
2375 call file_set_attribute( fid,
"global",
"scale_cartesC_prc_num_x", (/prc_num_x/) )
2376 call file_set_attribute( fid,
"global",
"scale_cartesC_prc_num_y", (/prc_num_y/) )
2378 call file_set_attribute( fid,
"global",
"scale_cartesC_prc_periodic_z", .false. )
2379 call file_set_attribute( fid,
"global",
"scale_cartesC_prc_periodic_x", prc_periodic_x )
2380 call file_set_attribute( fid,
"global",
"scale_cartesC_prc_periodic_y", prc_periodic_y )
2382 call file_set_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_imaxg", (/
imaxg/) )
2383 call file_set_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_jmaxg", (/
jmaxg/) )
2385 call file_set_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_kmax", (/
kmax/) )
2386 if (
okmax > 0 )
call file_set_attribute( fid,
"global",
"scale_ocean_grid_cartesC_index_kmax", (/
okmax/) )
2387 if (
lkmax > 0 )
call file_set_attribute( fid,
"global",
"scale_land_grid_cartesC_index_kmax", (/
lkmax/) )
2388 if (
ukmax > 0 )
call file_set_attribute( fid,
"global",
"scale_urban_grid_cartesC_index_kmax", (/
ukmax/) )
2390 call file_set_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_khalo", (/
khalo/) )
2391 call file_set_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_ihalo", (/
ihalo/) )
2392 call file_set_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_jhalo", (/
jhalo/) )
2394 if ( calendar /=
"" )
call file_set_attribute( fid,
"global",
"calendar", calendar )
2395 call file_set_attribute( fid,
"global",
"time_units", tunits )
2396 call file_set_attribute( fid,
"global",
"time_start", (/time/) )
2413 file_set_attribute, &
2422 integer,
intent(in) :: fid
2423 integer,
intent(in) :: dtype
2424 logical,
intent(in) :: hasZ
2431 character(len=2) :: axisname(3)
2433 logical,
save :: set_dim = .false.
2438 if ( .not. set_dim )
then
2457 call file_def_axis( fid,
'zh' ,
'Z (half level)' ,
'm',
'zh' , dtype,
kmax+1, bounds=.true. )
2459 if (
okmax > 0 )
then
2461 call file_def_axis( fid,
'ozh',
'OZ (half level)',
'm',
'ozh', dtype,
okmax+1, bounds=.true. )
2464 if (
lkmax > 0 )
then
2466 call file_def_axis( fid,
'lzh',
'LZ (half level)',
'm',
'lzh', dtype,
lkmax+1, bounds=.true. )
2469 if (
ukmax > 0 )
then
2471 call file_def_axis( fid,
'uzh',
'UZ (half level)',
'm',
'uzh', dtype,
ukmax+1, bounds=.true. )
2475 call file_def_axis( fid,
'x' ,
'X' ,
'm',
'x' , dtype, isize, bounds=.true. )
2476 call file_def_axis( fid,
'xh' ,
'X (half level)' ,
'm',
'xh' , dtype, isize, bounds=.true. )
2477 call file_def_axis( fid,
'y' ,
'Y' ,
'm',
'y' , dtype, jsize, bounds=.true. )
2478 call file_def_axis( fid,
'yh' ,
'Y (half level)' ,
'm',
'yh' , dtype, jsize, bounds=.true. )
2481 call file_def_axis( fid,
'CZ' ,
'Atmos Grid Center Position Z',
'm',
'CZ', dtype,
ka )
2482 call file_def_axis( fid,
'FZ' ,
'Atmos Grid Face Position Z',
'm',
'FZ', dtype,
ka+1 )
2483 call file_def_axis( fid,
'CDZ' ,
'Grid Cell length Z',
'm',
'CZ', dtype,
ka )
2484 call file_def_axis( fid,
'FDZ' ,
'Grid distance Z',
'm',
'FDZ', dtype,
ka-1 )
2485 call file_def_axis( fid,
'CBFZ' ,
'Boundary factor Center Z',
'1',
'CZ', dtype,
ka )
2486 call file_def_axis( fid,
'FBFZ' ,
'Boundary factor Face Z',
'1',
'FZ', dtype,
ka+1 )
2488 if (
okmax > 0 )
then
2489 call file_def_axis( fid,
'OCZ' ,
'Ocean Grid Center Position Z',
'm',
'OCZ', dtype,
okmax )
2490 call file_def_axis( fid,
'OFZ' ,
'Ocean Grid Face Position Z',
'm',
'OFZ', dtype,
okmax+1 )
2491 call file_def_axis( fid,
'OCDZ' ,
'Ocean Grid Cell length Z',
'm',
'OCZ', dtype,
okmax )
2494 if (
lkmax > 0 )
then
2495 call file_def_axis( fid,
'LCZ' ,
'Land Grid Center Position Z',
'm',
'LCZ', dtype,
lkmax )
2496 call file_def_axis( fid,
'LFZ' ,
'Land Grid Face Position Z',
'm',
'LFZ', dtype,
lkmax+1 )
2497 call file_def_axis( fid,
'LCDZ' ,
'Land Grid Cell length Z',
'm',
'LCZ', dtype,
lkmax )
2500 if (
ukmax > 0 )
then
2501 call file_def_axis( fid,
'UCZ' ,
'Urban Grid Center Position Z',
'm',
'UCZ', dtype,
ukmax )
2502 call file_def_axis( fid,
'UFZ' ,
'Urban Grid Face Position Z',
'm',
'UFZ', dtype,
ukmax+1 )
2503 call file_def_axis( fid,
'UCDZ' ,
'Urban Grid Cell length Z',
'm',
'UCZ', dtype,
ukmax )
2507 call file_def_axis( fid,
'CX' ,
'Atmos Grid Center Position X',
'm',
'CX', dtype, iall )
2508 call file_def_axis( fid,
'CY' ,
'Atmos Grid Center Position Y',
'm',
'CY', dtype, jall )
2509 call file_def_axis( fid,
'FX' ,
'Atmos Grid Face Position X',
'm',
'FX', dtype, iall+1 )
2510 call file_def_axis( fid,
'FY' ,
'Atmos Grid Face Position Y',
'm',
'FY', dtype, jall+1 )
2511 call file_def_axis( fid,
'CDX' ,
'Grid Cell length X',
'm',
'CX', dtype, iall )
2512 call file_def_axis( fid,
'CDY' ,
'Grid Cell length Y',
'm',
'CY', dtype, jall )
2513 call file_def_axis( fid,
'FDX' ,
'Grid distance X',
'm',
'FX', dtype, iall+1 )
2514 call file_def_axis( fid,
'FDY' ,
'Grid distance Y',
'm',
'FY', dtype, jall+1 )
2515 call file_def_axis( fid,
'CBFX' ,
'Boundary factor Center X',
'1',
'CX', dtype, iall )
2516 call file_def_axis( fid,
'CBFY' ,
'Boundary factor Center Y',
'1',
'CY', dtype, jall )
2517 call file_def_axis( fid,
'FBFX' ,
'Boundary factor Face X',
'1',
'FX', dtype, iall+1 )
2518 call file_def_axis( fid,
'FBFY' ,
'Boundary factor Face Y',
'1',
'FY', dtype, jall+1 )
2520 call file_def_axis( fid,
'CXG' ,
'Grid Center Position X (global)',
'm',
'CXG', dtype,
iag )
2521 call file_def_axis( fid,
'CYG' ,
'Grid Center Position Y (global)',
'm',
'CYG', dtype,
jag )
2522 call file_def_axis( fid,
'FXG' ,
'Grid Face Position X (global)',
'm',
'FXG', dtype,
iag+1 )
2523 call file_def_axis( fid,
'FYG' ,
'Grid Face Position Y (global)',
'm',
'FYG', dtype,
jag+1 )
2524 call file_def_axis( fid,
'CDXG' ,
'Grid Cell length X (global)',
'm',
'CXG', dtype,
iag )
2525 call file_def_axis( fid,
'CDYG' ,
'Grid Cell length Y (global)',
'm',
'CYG', dtype,
jag )
2526 call file_def_axis( fid,
'FDXG' ,
'Grid distance X (global)',
'm',
'FXG', dtype,
iag+1 )
2527 call file_def_axis( fid,
'FDYG' ,
'Grid distance Y (global)',
'm',
'FYG', dtype,
jag+1 )
2528 call file_def_axis( fid,
'CBFXG',
'Boundary factor Center X (global)',
'1',
'CXG', dtype,
iag )
2529 call file_def_axis( fid,
'CBFYG',
'Boundary factor Center Y (global)',
'1',
'CYG', dtype,
jag )
2530 call file_def_axis( fid,
'FBFXG',
'Boundary factor Face X (global)',
'1',
'FXG', dtype,
iag+1 )
2531 call file_def_axis( fid,
'FBFYG',
'Boundary factor Face Y (global)',
'1',
'FYG', dtype,
jag+1 )
2534 axisname(1:2) = (/
'x ',
'y '/)
2536 axisname(1:2) = (/
'xh',
'y '/)
2538 axisname(1:2) = (/
'x ',
'yh'/)
2540 axisname(1:2) = (/
'xh',
'yh'/)
2542 axisname(1:2) = (/
'x ',
'y '/)
2544 axisname(1:2) = (/
'xh',
'y '/)
2546 axisname(1:2) = (/
'x ',
'yh'/)
2548 axisname(1:2) = (/
'xh',
'yh'/)
2552 axisname(1:2) = (/
'x ',
'y '/)
2555 axisname(1:2) = (/
'x ',
'y '/)
2558 axisname(1:2) = (/
'x ',
'y '/)
2560 axisname(1:2) = (/
'xh',
'y '/)
2562 axisname(1:2) = (/
'x ',
'yh'/)
2564 axisname(1:2) = (/
'xh',
'yh'/)
2567 axisname = (/
'z ',
'x ',
'y '/)
2569 'm', axisname(1:3), dtype )
2570 axisname = (/
'zh',
'x ',
'y '/)
2572 'm', axisname(1:3), dtype )
2574 axisname = (/
'z ',
'xh',
'y '/)
2576 'm2', axisname(1:3), dtype )
2577 axisname = (/
'z ',
'x ',
'yh'/)
2579 'm2', axisname(1:3), dtype )
2580 axisname = (/
'zh',
'xh',
'y '/)
2582 'm2', axisname(1:3), dtype )
2583 axisname = (/
'zh',
'x ',
'yh'/)
2585 'm2', axisname(1:3), dtype )
2586 axisname = (/
'z ',
'x ',
'y '/)
2588 'm2', axisname(1:3), dtype )
2589 axisname = (/
'z ',
'xh',
'yh'/)
2591 'm2', axisname(1:3), dtype )
2592 axisname = (/
'z ',
'xh',
'yh'/)
2594 'm2', axisname(1:3), dtype )
2595 axisname = (/
'z ',
'x ',
'y '/)
2597 'm2', axisname(1:3), dtype )
2599 axisname = (/
'z ',
'x ',
'y '/)
2601 'm3', axisname(1:3), dtype )
2602 axisname = (/
'zh',
'x ',
'y '/)
2604 'm3', axisname(1:3), dtype )
2605 axisname = (/
'z ',
'xh',
'y '/)
2607 'm3', axisname(1:3), dtype )
2608 axisname = (/
'z ',
'x ',
'yh'/)
2610 'm3', axisname(1:3), dtype )
2612 if (
okmax > 0 )
then
2613 axisname = (/
'oz',
'x ',
'y '/)
2615 'm3', axisname(1:3), dtype )
2617 if (
lkmax > 0 )
then
2618 axisname = (/
'lz',
'x ',
'y '/)
2620 'm3', axisname(1:3), dtype )
2622 if (
ukmax > 0 )
then
2623 axisname = (/
'uz',
'x ',
'y '/)
2625 'm3', axisname(1:3), dtype )
2632 if (
okmax > 0 )
then
2633 call file_set_attribute( fid,
'oz' ,
'positive',
'down' )
2634 call file_set_attribute( fid,
'ozh',
'positive',
'down' )
2636 if (
lkmax > 0 )
then
2637 call file_set_attribute( fid,
'lz' ,
'positive',
'down' )
2638 call file_set_attribute( fid,
'lzh',
'positive',
'down' )
2640 if (
ukmax > 0 )
then
2641 call file_set_attribute( fid,
'uz' ,
'positive',
'down' )
2642 call file_set_attribute( fid,
'uzh',
'positive',
'down' )
2644 if (
okmax > 0 )
then
2645 call file_set_attribute( fid,
'OCZ',
'positive',
'down' )
2646 call file_set_attribute( fid,
'OFZ',
'positive',
'down' )
2648 if (
lkmax > 0 )
then
2649 call file_set_attribute( fid,
'LCZ',
'positive',
'down' )
2650 call file_set_attribute( fid,
'LFZ',
'positive',
'down' )
2652 if (
ukmax > 0 )
then
2653 call file_set_attribute( fid,
'UCZ',
'positive',
'down' )
2654 call file_set_attribute( fid,
'UFZ',
'positive',
'down' )
2659 call file_set_attribute( fid,
"x" ,
"size_global" , (/
iag /) )
2660 call file_set_attribute( fid,
"x" ,
"start_global", (/ 1 /) )
2661 call file_set_attribute( fid,
"x" ,
"halo_global" , (/
ihalo,
ihalo /) )
2662 call file_set_attribute( fid,
"x" ,
"halo_local" , (/
ihalo,
ihalo /) )
2664 call file_set_attribute( fid,
"xh",
"size_global" , (/
iag+1 /) )
2665 call file_set_attribute( fid,
"xh",
"start_global", (/ 1 /) )
2666 call file_set_attribute( fid,
"xh",
"halo_global" , (/
ihalo,
ihalo /) )
2667 call file_set_attribute( fid,
"xh",
"halo_local" , (/
ihalo,
ihalo /) )
2669 call file_set_attribute( fid,
"y" ,
"size_global" , (/
jag /) )
2670 call file_set_attribute( fid,
"y" ,
"start_global", (/ 1 /) )
2671 call file_set_attribute( fid,
"y" ,
"halo_global" , (/
jhalo,
jhalo /) )
2672 call file_set_attribute( fid,
"y" ,
"halo_local" , (/
jhalo,
jhalo /) )
2674 call file_set_attribute( fid,
"yh",
"size_global" , (/
jag+1 /) )
2675 call file_set_attribute( fid,
"yh",
"start_global", (/ 1 /) )
2676 call file_set_attribute( fid,
"yh",
"halo_global" , (/
jhalo,
jhalo /) )
2677 call file_set_attribute( fid,
"yh",
"halo_local" , (/
jhalo,
jhalo /) )
2679 call file_set_attribute( fid,
"x" ,
"size_global" , file_cartesc_axis_info(1)%size_global (:) )
2680 call file_set_attribute( fid,
"x" ,
"start_global", file_cartesc_axis_info(1)%start_global(:) )
2681 call file_set_attribute( fid,
"x" ,
"halo_global" , file_cartesc_axis_info(1)%halo_global (:) )
2682 call file_set_attribute( fid,
"x" ,
"halo_local" , file_cartesc_axis_info(1)%halo_local (:) )
2684 call file_set_attribute( fid,
"xh",
"size_global" , file_cartesc_axis_info(2)%size_global (:) )
2685 call file_set_attribute( fid,
"xh",
"start_global", file_cartesc_axis_info(2)%start_global(:) )
2686 call file_set_attribute( fid,
"xh",
"halo_global" , file_cartesc_axis_info(2)%halo_global (:) )
2687 call file_set_attribute( fid,
"xh",
"halo_local" , file_cartesc_axis_info(2)%halo_local (:) )
2689 call file_set_attribute( fid,
"y" ,
"size_global" , file_cartesc_axis_info(3)%size_global (:) )
2690 call file_set_attribute( fid,
"y" ,
"start_global", file_cartesc_axis_info(3)%start_global(:) )
2691 call file_set_attribute( fid,
"y" ,
"halo_global" , file_cartesc_axis_info(3)%halo_global (:) )
2692 call file_set_attribute( fid,
"y" ,
"halo_local" , file_cartesc_axis_info(3)%halo_local (:) )
2694 call file_set_attribute( fid,
"yh",
"size_global" , file_cartesc_axis_info(4)%size_global (:) )
2695 call file_set_attribute( fid,
"yh",
"start_global", file_cartesc_axis_info(4)%start_global(:) )
2696 call file_set_attribute( fid,
"yh",
"halo_global" , file_cartesc_axis_info(4)%halo_global (:) )
2697 call file_set_attribute( fid,
"yh",
"halo_local" , file_cartesc_axis_info(4)%halo_local (:) )
2700 call file_set_attribute( fid,
"x" ,
"periodic" , file_cartesc_axis_info(1)%periodic )
2701 call file_set_attribute( fid,
"xh",
"periodic" , file_cartesc_axis_info(2)%periodic )
2702 call file_set_attribute( fid,
"y" ,
"periodic" , file_cartesc_axis_info(3)%periodic )
2703 call file_set_attribute( fid,
"yh",
"periodic" , file_cartesc_axis_info(4)%periodic )
2708 call file_set_attribute( fid,
"x" ,
"standard_name",
"projection_x_coordinate" )
2709 call file_set_attribute( fid,
"xh",
"standard_name",
"projection_x_coordinate" )
2710 call file_set_attribute( fid,
"y" ,
"standard_name",
"projection_y_coordinate" )
2711 call file_set_attribute( fid,
"yh",
"standard_name",
"projection_y_coordinate" )
2717 call file_set_attribute( fid, &
2724 call file_set_attribute( fid, &
2731 call file_set_attribute( fid, &
2733 "longitude_of_central_meridian", &
2738 call file_set_attribute( fid, &
2740 "longitude_of_projection_origin", &
2745 call file_set_attribute( fid, &
2747 "latitude_of_projection_origin", &
2752 call file_set_attribute( fid, &
2754 "straight_vertical_longitude_from_pole", &
2760 call file_set_attribute( fid, &
2762 "standard_parallel", &
2765 call file_set_attribute( fid, &
2767 "standard_parallel", &
2773 call file_set_attribute( fid, &
2783 call file_set_attribute( fid,
"cell_area",
"standard_name",
"area" )
2784 call file_set_attribute( fid,
"cell_area_uy",
"standard_name",
"area" )
2785 call file_set_attribute( fid,
"cell_area_xv",
"standard_name",
"area" )
2788 call file_set_attribute( fid,
"cell_area_zuy_x",
"standard_name",
"area" )
2789 call file_set_attribute( fid,
"cell_area_zxv_y",
"standard_name",
"area" )
2790 call file_set_attribute( fid,
"cell_area_wuy_x",
"standard_name",
"area" )
2791 call file_set_attribute( fid,
"cell_area_wxv_y",
"standard_name",
"area" )
2792 call file_set_attribute( fid,
"cell_area_zxy_x",
"standard_name",
"area" )
2793 call file_set_attribute( fid,
"cell_area_zuv_y",
"standard_name",
"area" )
2794 call file_set_attribute( fid,
"cell_area_zuv_x",
"standard_name",
"area" )
2795 call file_set_attribute( fid,
"cell_area_zxy_y",
"standard_name",
"area" )
2797 call file_set_attribute( fid,
"cell_volume",
"standard_name",
"volume" )
2798 call file_set_attribute( fid,
"cell_volume_wxy",
"standard_name",
"volume" )
2799 call file_set_attribute( fid,
"cell_volume_zuy",
"standard_name",
"volume" )
2800 call file_set_attribute( fid,
"cell_volume_zxv",
"standard_name",
"volume" )
2802 if (
okmax > 0 )
then
2803 call file_set_attribute( fid,
"cell_volume_oxy",
"standard_name",
"volume" )
2805 if (
lkmax > 0 )
then
2806 call file_set_attribute( fid,
"cell_volume_lxy",
"standard_name",
"volume" )
2808 if (
ukmax > 0 )
then
2809 call file_set_attribute( fid,
"cell_volume_uxy",
"standard_name",
"volume" )
2815 call file_set_attribute( fid,
"grid",
"cf_role",
"grid_topology" )
2816 call file_set_attribute( fid,
"grid",
"topology_dimension", (/ 2 /) )
2817 call file_set_attribute( fid,
"grid",
"node_dimensions",
"xh yh" )
2818 call file_set_attribute( fid,
"grid",
"face_dimensions",
"x: xh (padding: none) y: yh (padding: none)" )
2819 call file_set_attribute( fid,
"grid",
"node_coordinates",
"lon_uv lat_uv" )
2820 call file_set_attribute( fid,
"grid",
"face_coordinates",
"lon lat" )
2821 call file_set_attribute( fid,
"grid",
"edge1_coordinates",
"lon_uy lat_uy" )
2822 call file_set_attribute( fid,
"grid",
"edge2_coordinates",
"lon_xv lat_xv" )
2823 call file_set_attribute( fid,
"grid",
"vertical_dimensions",
"z: zh (padding: none)" )
2825 if (
okmax > 0 )
then
2827 call file_set_attribute( fid,
"grid_ocean",
"cf_role",
"grid_topology" )
2828 call file_set_attribute( fid,
"grid_ocean",
"topology_dimension", (/ 2 /) )
2829 call file_set_attribute( fid,
"grid_ocean",
"node_dimensions",
"xh yh" )
2830 call file_set_attribute( fid,
"grid_ocean",
"face_dimensions",
"x: xh (padding: none) y: yh (padding: none)" )
2831 call file_set_attribute( fid,
"grid_ocean",
"node_coordinates",
"lon_uv lat_uv" )
2832 call file_set_attribute( fid,
"grid_ocean",
"face_coordinates",
"lon lat" )
2833 call file_set_attribute( fid,
"grid_ocean",
"edge1_coordinates",
"lon_uy lat_uy" )
2834 call file_set_attribute( fid,
"grid_ocean",
"edge2_coordinates",
"lon_xv lat_xv" )
2835 call file_set_attribute( fid,
"grid_ocean",
"vertical_dimensions",
"oz: ozh (padding: none)" )
2838 if (
lkmax > 0 )
then
2840 call file_set_attribute( fid,
"grid_land",
"cf_role",
"grid_topology" )
2841 call file_set_attribute( fid,
"grid_land",
"topology_dimension", (/ 2 /) )
2842 call file_set_attribute( fid,
"grid_land",
"node_dimensions",
"xh yh" )
2843 call file_set_attribute( fid,
"grid_land",
"face_dimensions",
"x: xh (padding: none) y: yh (padding: none)" )
2844 call file_set_attribute( fid,
"grid_land",
"node_coordinates",
"lon_uv lat_uv" )
2845 call file_set_attribute( fid,
"grid_land",
"face_coordinates",
"lon lat" )
2846 call file_set_attribute( fid,
"grid_land",
"edge1_coordinates",
"lon_uy lat_uy" )
2847 call file_set_attribute( fid,
"grid_land",
"edge2_coordinates",
"lon_xv lat_xv" )
2848 call file_set_attribute( fid,
"grid_land",
"vertical_dimensions",
"lz: lzh (padding: none)" )
2851 if (
ukmax > 0 )
then
2853 call file_set_attribute( fid,
"grid_urban",
"cf_role",
"grid_topology" )
2854 call file_set_attribute( fid,
"grid_urban",
"topology_dimension", (/ 2 /) )
2855 call file_set_attribute( fid,
"grid_urban",
"node_dimensions",
"xh yh" )
2856 call file_set_attribute( fid,
"grid_urban",
"face_dimensions",
"x: xh (padding: none) y: yh (padding: none)" )
2857 call file_set_attribute( fid,
"grid_urban",
"node_coordinates",
"lon_uv lat_uv" )
2858 call file_set_attribute( fid,
"grid_urban",
"face_coordinates",
"lon lat" )
2859 call file_set_attribute( fid,
"grid_urban",
"edge1_coordinates",
"lon_uy lat_uy" )
2860 call file_set_attribute( fid,
"grid_urban",
"edge2_coordinates",
"lon_xv lat_xv" )
2861 call file_set_attribute( fid,
"grid_urban",
"vertical_dimensions",
"uz: uzh (padding: none)" )
2865 call file_set_attribute( fid,
"grid_model",
"cf_role",
"grid_topology" )
2866 call file_set_attribute( fid,
"grid_model",
"topology_dimension", (/ 2 /) )
2867 call file_set_attribute( fid,
"grid_model",
"node_dimensions",
"FX FY" )
2868 call file_set_attribute( fid,
"grid_model",
"face_dimensions",
"CX: FY (padding: none) CY: FY (padding: none)" )
2869 call file_set_attribute( fid,
"grid_model",
"vertical_dimensions",
"CZ: FZ (padding: none)" )
2872 call file_set_attribute( fid,
"grid_model_global",
"cf_role",
"grid_topology" )
2873 call file_set_attribute( fid,
"grid_model_global",
"topology_dimension", (/ 2 /) )
2874 call file_set_attribute( fid,
"grid_model_global",
"node_dimensions",
"FXG FYG" )
2875 call file_set_attribute( fid,
"grid_model_global",
"face_dimensions",
"CXG: FYG (padding: none) CYG: FYG (padding: none)" )
2876 call file_set_attribute( fid,
"grid_model_global",
"vertical_dimensions",
"CZ: FZ (padding: none)" )
2893 file_write_associatedcoordinate
2944 integer,
intent(in) :: fid
2945 logical,
intent(in) :: haszcoord
2946 integer,
intent(in) :: start(3)
2948 logical :: put_z, put_x, put_y
2949 integer :: XSB, XEB, YSB, YEB
2951 real(RP) :: z_bnds(2,KA), zh_bnds(2,0:KA)
2952 real(RP) :: oz_bnds(2,OKA), ozh_bnds(2,0:OKA)
2953 real(RP) :: lz_bnds(2,LKA), lzh_bnds(2,0:LKA)
2954 real(RP) :: uz_bnds(2,UKA), uzh_bnds(2,0:UKA)
2955 real(RP) :: x_bnds(2,IA), xh_bnds(2,0:IA)
2956 real(RP) :: y_bnds(2,JA), yh_bnds(2,0:JA)
2957 integer :: start_(2)
2959 real(RP) :: FDXG(0:IAG), FDYG(0:JAG)
2960 real(RP) :: FDX(0:IA), FDY(0:JA)
2988 if ( haszcoord .and. put_z )
then
2990 start_(2) = start(1)
2998 call file_write_associatedcoordinate( fid,
'z_bnds', z_bnds(:,
ks:
ke), start_(:) )
3005 call file_write_associatedcoordinate( fid,
'zh_bnds', zh_bnds(:,
ks-1:
ke), start_(:) )
3008 if (
okmax > 0 )
then
3014 call file_write_associatedcoordinate( fid,
'oz_bnds', oz_bnds(:,
oks:
oke), start_(:) )
3023 call file_write_associatedcoordinate( fid,
'ozh_bnds', ozh_bnds(:,
oks-1:
oke), start_(:) )
3027 if (
lkmax > 0 )
then
3033 call file_write_associatedcoordinate( fid,
'lz_bnds', lz_bnds(:,
lks:
lke), start_(:) )
3042 call file_write_associatedcoordinate( fid,
'lzh_bnds', lzh_bnds(:,
lks-1:
lke), start_(:) )
3046 if (
ukmax > 0 )
then
3052 call file_write_associatedcoordinate( fid,
'uz_bnds', uz_bnds(:,
uks:
uke), start_(:) )
3061 call file_write_associatedcoordinate( fid,
'uzh_bnds', uzh_bnds(:,
uks-1:
uke), start_(:) )
3067 start_(2) = start(2)
3075 call file_write_associatedcoordinate( fid,
'x_bnds', x_bnds(:,isb2:ieb2), start_(:) )
3083 if ( ieb2 == ia )
then
3088 call file_write_associatedcoordinate( fid,
'xh_bnds', xh_bnds(:,isb2:ieb2), start_(:) )
3095 call file_write_associatedcoordinate( fid,
'x_bnds', x_bnds(:,
isb:
ieb), start_(:) )
3103 if (
ieb == ia )
then
3108 call file_write_associatedcoordinate( fid,
'xh_bnds', xh_bnds(:,
isb:
ieb), start_(:) )
3114 start_(2) = start(3)
3122 call file_write_associatedcoordinate( fid,
'y_bnds', y_bnds(:,jsb2:jeb2), start_(:) )
3130 if ( jeb2 == ja )
then
3135 call file_write_associatedcoordinate( fid,
'yh_bnds', yh_bnds(:,jsb2:jeb2), start_(:) )
3142 call file_write_associatedcoordinate( fid,
'y_bnds', y_bnds(:,
jsb:
jeb), start_(:) )
3150 if (
jeb == ja )
then
3155 call file_write_associatedcoordinate( fid,
'yh_bnds', yh_bnds(:,
jsb:
jeb), start_(:) )
3174 if ( haszcoord .and. put_z )
then
3182 if (
okmax > 0 )
then
3188 if (
lkmax > 0 )
then
3194 if (
ukmax > 0 )
then
3210 call file_write_axis( fid,
'FDX', fdxg(:) )
3211 call file_write_axis( fid,
'FDY', fdyg(:) )
3224 call file_write_axis( fid,
'FDX', fdx(:) )
3225 call file_write_axis( fid,
'FDY', fdy(:) )
3238 call file_write_axis( fid,
'FDXG', fdxg(:) )
3239 call file_write_axis( fid,
'FDYG', fdyg(:) )
3247 call file_write_associatedcoordinate( fid,
'lon' , axis_lon(:,:), start(2:3) )
3248 call file_write_associatedcoordinate( fid,
'lon_uy', axis_lonuy(:,:), start(2:3) )
3249 call file_write_associatedcoordinate( fid,
'lon_xv', axis_lonxv(:,:), start(2:3) )
3250 call file_write_associatedcoordinate( fid,
'lon_uv', axis_lonuv(:,:), start(2:3) )
3251 call file_write_associatedcoordinate( fid,
'lat' , axis_lat(:,:), start(2:3) )
3252 call file_write_associatedcoordinate( fid,
'lat_uy', axis_latuy(:,:), start(2:3) )
3253 call file_write_associatedcoordinate( fid,
'lat_xv', axis_latxv(:,:), start(2:3) )
3254 call file_write_associatedcoordinate( fid,
'lat_uv', axis_latuv(:,:), start(2:3) )
3256 if ( haszcoord )
then
3257 call file_write_associatedcoordinate( fid,
'topo', axis_topo(:,:), start(2:3) )
3259 call file_write_associatedcoordinate( fid,
'lsmask', axis_lsmask(:,:), start(2:3) )
3261 call file_write_associatedcoordinate( fid,
'cell_area', axis_area(:,:), start(2:3) )
3262 call file_write_associatedcoordinate( fid,
'cell_area_uy', axis_areauy(:,:), start(2:3) )
3263 call file_write_associatedcoordinate( fid,
'cell_area_xv', axis_areaxv(:,:), start(2:3) )
3265 if ( haszcoord )
then
3266 call file_write_associatedcoordinate( fid,
'height' , axis_hgt(:,:,:), start(1:3) )
3267 call file_write_associatedcoordinate( fid,
'height_wxy', axis_hgtwxy(:,:,:), start(1:3) )
3269 call file_write_associatedcoordinate( fid,
'cell_area_zuy_x', axis_areazuy_x(:,:,:), start(1:3) )
3270 call file_write_associatedcoordinate( fid,
'cell_area_zxv_y', axis_areazxv_y(:,:,:), start(1:3) )
3271 call file_write_associatedcoordinate( fid,
'cell_area_wuy_x', axis_areawuy_x(:,:,:), start(1:3) )
3272 call file_write_associatedcoordinate( fid,
'cell_area_wxv_y', axis_areawxv_y(:,:,:), start(1:3) )
3273 call file_write_associatedcoordinate( fid,
'cell_area_zxy_x', axis_areazxy_x(:,:,:), start(1:3) )
3274 call file_write_associatedcoordinate( fid,
'cell_area_zuv_y', axis_areazuv_y(:,:,:), start(1:3) )
3275 call file_write_associatedcoordinate( fid,
'cell_area_zuv_x', axis_areazuv_x(:,:,:), start(1:3) )
3276 call file_write_associatedcoordinate( fid,
'cell_area_zxy_y', axis_areazxy_y(:,:,:), start(1:3) )
3278 call file_write_associatedcoordinate( fid,
'cell_volume', axis_vol(:,:,:), start(1:3) )
3279 call file_write_associatedcoordinate( fid,
'cell_volume_wxy', axis_volwxy(:,:,:), start(1:3) )
3280 call file_write_associatedcoordinate( fid,
'cell_volume_zuy', axis_volzuy(:,:,:), start(1:3) )
3281 call file_write_associatedcoordinate( fid,
'cell_volume_zxv', axis_volzxv(:,:,:), start(1:3) )
3283 if (
okmax > 0 )
then
3284 call file_write_associatedcoordinate( fid,
'cell_volume_oxy', axis_volo(:,:,:), start(1:3) )
3286 if (
lkmax > 0 )
then
3287 call file_write_associatedcoordinate( fid,
'cell_volume_lxy', axis_voll(:,:,:), start(1:3) )
3289 if (
ukmax > 0 )
then
3290 call file_write_associatedcoordinate( fid,
'cell_volume_uxy', axis_volu(:,:,:), start(1:3) )
3294 xsb =
isb - isb2 + 1
3296 ysb =
jsb - jsb2 + 1
3299 call file_write_associatedcoordinate( fid,
'lon' , axis_lon(xsb:xeb,ysb:yeb), start(2:3) )
3300 call file_write_associatedcoordinate( fid,
'lon_uy', axis_lonuy(xsb:xeb,ysb:yeb), start(2:3) )
3301 call file_write_associatedcoordinate( fid,
'lon_xv', axis_lonxv(xsb:xeb,ysb:yeb), start(2:3) )
3302 call file_write_associatedcoordinate( fid,
'lon_uv', axis_lonuv(xsb:xeb,ysb:yeb), start(2:3) )
3303 call file_write_associatedcoordinate( fid,
'lat' , axis_lat(xsb:xeb,ysb:yeb), start(2:3) )
3304 call file_write_associatedcoordinate( fid,
'lat_uy', axis_latuy(xsb:xeb,ysb:yeb), start(2:3) )
3305 call file_write_associatedcoordinate( fid,
'lat_xv', axis_latxv(xsb:xeb,ysb:yeb), start(2:3) )
3306 call file_write_associatedcoordinate( fid,
'lat_uv', axis_latuv(xsb:xeb,ysb:yeb), start(2:3) )
3308 if ( haszcoord )
then
3309 call file_write_associatedcoordinate( fid,
'topo', axis_topo(xsb:xeb,ysb:yeb), start(2:3) )
3311 call file_write_associatedcoordinate( fid,
'lsmask', axis_lsmask(xsb:xeb,ysb:yeb), start(2:3) )
3313 call file_write_associatedcoordinate( fid,
'cell_area', axis_area(xsb:xeb,ysb:yeb), start(2:3) )
3314 call file_write_associatedcoordinate( fid,
'cell_area_uy', axis_areauy(xsb:xeb,ysb:yeb), start(2:3) )
3315 call file_write_associatedcoordinate( fid,
'cell_area_xv', axis_areaxv(xsb:xeb,ysb:yeb), start(2:3) )
3317 if ( haszcoord )
then
3318 call file_write_associatedcoordinate( fid,
'height' , axis_hgt(:,xsb:xeb,ysb:yeb), start(1:3) )
3319 call file_write_associatedcoordinate( fid,
'height_wxy', axis_hgtwxy(:,xsb:xeb,ysb:yeb), start(1:3) )
3321 call file_write_associatedcoordinate( fid,
'cell_area_zuy_x', axis_areazuy_x(:,xsb:xeb,ysb:yeb), start(1:3) )
3322 call file_write_associatedcoordinate( fid,
'cell_area_zxv_y', axis_areazxv_y(:,xsb:xeb,ysb:yeb), start(1:3) )
3323 call file_write_associatedcoordinate( fid,
'cell_area_wuy_x', axis_areawuy_x(:,xsb:xeb,ysb:yeb), start(1:3) )
3324 call file_write_associatedcoordinate( fid,
'cell_area_wxv_y', axis_areawxv_y(:,xsb:xeb,ysb:yeb), start(1:3) )
3325 call file_write_associatedcoordinate( fid,
'cell_area_zxy_x', axis_areazxy_x(:,xsb:xeb,ysb:yeb), start(1:3) )
3326 call file_write_associatedcoordinate( fid,
'cell_area_zuv_y', axis_areazuv_y(:,xsb:xeb,ysb:yeb), start(1:3) )
3327 call file_write_associatedcoordinate( fid,
'cell_area_zuv_x', axis_areazuv_x(:,xsb:xeb,ysb:yeb), start(1:3) )
3328 call file_write_associatedcoordinate( fid,
'cell_area_zxy_y', axis_areazxy_y(:,xsb:xeb,ysb:yeb), start(1:3) )
3330 call file_write_associatedcoordinate( fid,
'cell_volume', axis_vol(:,xsb:xeb,ysb:yeb), start(1:3) )
3331 call file_write_associatedcoordinate( fid,
'cell_volume_wxy', axis_volwxy(:,xsb:xeb,ysb:yeb), start(1:3) )
3332 call file_write_associatedcoordinate( fid,
'cell_volume_zuy', axis_volzuy(:,xsb:xeb,ysb:yeb), start(1:3) )
3333 call file_write_associatedcoordinate( fid,
'cell_volume_zxv', axis_volzxv(:,xsb:xeb,ysb:yeb), start(1:3) )
3335 if (
okmax > 0 )
then
3336 call file_write_associatedcoordinate( fid,
'cell_volume_oxy', axis_volo(:,xsb:xeb,ysb:yeb), start(1:3) )
3338 if (
lkmax > 0 )
then
3339 call file_write_associatedcoordinate( fid,
'cell_volume_lxy', axis_voll(:,xsb:xeb,ysb:yeb), start(1:3) )
3341 if (
ukmax > 0 )
then
3342 call file_write_associatedcoordinate( fid,
'cell_volume_uxy', axis_volu(:,xsb:xeb,ysb:yeb), start(1:3) )
3354 varname, desc, unit, &
3355 dim_type, datatype, &
3376 integer,
intent(in) :: fid
3377 character(len=*),
intent(in) :: varname
3378 character(len=*),
intent(in) :: desc
3379 character(len=*),
intent(in) :: unit
3380 character(len=*),
intent(in) :: dim_type
3381 character(len=*),
intent(in) :: datatype
3383 integer,
intent(out) :: vid
3385 character(len=*),
intent(in),
optional :: standard_name
3386 real(
dp),
intent(in),
optional :: timeintv
3387 integer,
intent(in),
optional :: nsteps
3388 character(len=*),
intent(in),
optional :: cell_measures
3390 character(len=H_MID) :: standard_name_
3391 character(len=H_SHORT) :: cell_measures_
3393 character(len=H_SHORT) :: dimtype
3395 integer :: dtype, elm_size, ndims
3403 if ( datatype ==
'REAL8' )
then
3406 elseif( datatype ==
'REAL4' )
then
3413 elseif(
rp == 4 )
then
3417 log_error(
"FILE_CARTESC_def_var",*)
'unsupported data type. Check!', trim(datatype),
' item:',trim(varname)
3423 select case( dim_type )
3442 do n = 1, file_cartesc_ndims
3443 if ( file_cartesc_dims(n)%name == dimtype )
then
3448 if ( dimid <= -1 )
then
3449 log_error(
"FILE_CARTESC_def_var",*)
'dim_type is not supported: ', trim(dimtype)
3453 if (
present(nsteps) )
then
3454 write_buf_amount(fid) = write_buf_amount(fid) + file_cartesc_dims(dimid)%size * elm_size * nsteps
3456 write_buf_amount(fid) = write_buf_amount(fid) + file_cartesc_dims(dimid)%size * elm_size
3459 ndims = file_cartesc_dims(dimid)%ndims
3461 if (
present(standard_name) )
then
3462 standard_name_ = standard_name
3467 if (
present(timeintv) )
then
3469 ndims, file_cartesc_dims(dimid)%dims(1:ndims), dtype, &
3474 ndims, file_cartesc_dims(dimid)%dims(1:ndims), dtype, &
3478 if (
present(cell_measures) )
then
3479 cell_measures_ = cell_measures
3480 select case ( cell_measures )
3482 if ( file_cartesc_dims(dimid)%area ==
"" )
then
3483 log_error(
"FILE_CARTESC_def_var",*)
'area is not supported for ', trim(dimtype),
' as cell_measures'
3487 if ( file_cartesc_dims(dimid)%area ==
"" )
then
3488 log_error(
"FILE_CARTESC_def_var",*)
'area_z is not supported for ', trim(dimtype),
' as cell_measures'
3492 if ( file_cartesc_dims(dimid)%area_x ==
"" )
then
3493 log_error(
"FILE_CARTESC_def_var",*)
'area_x is not supported for ', trim(dimtype),
' as cell_measures'
3497 if ( file_cartesc_dims(dimid)%area_y ==
"" )
then
3498 log_error(
"FILE_CARTESC_def_var",*)
'area_y is not supported for ', trim(dimtype),
' as cell_measures'
3502 if ( file_cartesc_dims(dimid)%volume ==
"" )
then
3503 log_error(
"FILE_CARTESC_def_var",*)
'volume is not supported for ', trim(dimtype),
' as cell_measures'
3507 log_error(
"FILE_CARTESC_def_var",*)
'cell_measures must be "area" or "volume"'
3510 else if ( ndims == 2 )
then
3511 cell_measures_ =
"area"
3512 else if ( ndims == 3 )
then
3513 cell_measures_ =
"volume"
3518 select case( cell_measures_ )
3519 case (
"area",
"area_z" )
3520 call file_set_attribute( fid, varname,
"cell_measures",
"area: "//trim(file_cartesc_dims(dimid)%area) )
3522 call file_set_attribute( fid, varname,
"cell_measures",
"area: "//trim(file_cartesc_dims(dimid)%area_x) )
3524 call file_set_attribute( fid, varname,
"cell_measures",
"area: "//trim(file_cartesc_dims(dimid)%area_y) )
3526 call file_set_attribute( fid, varname,
"cell_measures",
"volume: "//trim(file_cartesc_dims(dimid)%volume) )
3535 if ( file_cartesc_dims(dimid)%location /=
"" )
then
3536 call file_set_attribute( fid, varname,
"grid", file_cartesc_dims(dimid)%grid )
3537 call file_set_attribute( fid, varname,
"location", file_cartesc_dims(dimid)%location )
3567 integer,
intent(in) :: fid
3568 integer,
intent(in) :: vid
3569 real(RP),
intent(in) :: var(:)
3570 character(len=*),
intent(in) :: varname
3571 character(len=*),
intent(in) :: dim_type
3573 integer :: dim1_S, dim1_E
3574 integer :: rankidx(2)
3586 if ( dim_type ==
'Z' )
then
3595 elseif( dim_type ==
'X' )
then
3597 exec = ( rankidx(2) == 0 )
3606 elseif( dim_type ==
'Y' )
then
3608 exec = ( rankidx(1) == 0 )
3618 log_error(
"FILE_CARTESC_write_var_1D",*)
'unsupported dimenstion type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
3622 if( exec )
call file_write( vid, var(dim1_s:dim1_e), &
3623 nowdaysec, nowdaysec, start=start )
3657 integer,
intent(in) :: fid
3658 integer,
intent(in) :: vid
3659 real(RP),
intent(in) :: var(:,:)
3660 character(len=*),
intent(in) :: varname
3661 character(len=*),
intent(in) :: dim_type
3662 logical,
intent(in),
optional :: fill_halo
3664 real(RP),
allocatable :: buf(:,:)
3666 integer :: dim1_S, dim1_E
3667 integer :: dim2_S, dim2_E
3670 logical :: fill_halo_
3671 integer :: rankidx(2)
3686 fill_halo_ = .false.
3687 if(
present(fill_halo) ) fill_halo_ = fill_halo
3689 if ( dim_type ==
'XY' &
3690 .OR. dim_type ==
'UY' &
3691 .OR. dim_type ==
'XV' &
3692 .OR. dim_type ==
'UV' )
then
3705 elseif( dim_type ==
'ZX' )
then
3711 exec = ( rankidx(2) == 0 )
3720 log_error(
"FILE_CARTESC_write_var_2D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
3728 allocate( buf(dim1_s:dim1_e,dim2_s:dim2_e) )
3730 if ( fill_halo_ )
then
3738 do j = dim2_s, dim2_e
3744 do j = dim2_s, dim2_e
3751 do i = dim1_s, dim1_e
3757 do i = dim1_s, dim1_e
3762 do j = dim2_s, dim2_e
3763 do i = dim1_s, dim1_e
3769 call file_write( vid, buf(:,:), nowdaysec, nowdaysec, start )
3808 integer,
intent(in) :: fid
3809 integer,
intent(in) :: vid
3810 real(RP),
intent(in) :: var(:,:,:)
3811 character(len=*),
intent(in) :: varname
3812 character(len=*),
intent(in) :: dim_type
3814 logical,
intent(in),
optional :: fill_halo
3816 real(RP),
allocatable :: buf(:,:,:)
3818 integer :: dim1_S, dim1_E
3819 integer :: dim2_S, dim2_E
3820 integer :: dim3_S, dim3_E
3823 logical :: fill_halo_
3824 integer :: rankidx(2)
3832 fill_halo_ = .false.
3833 if(
present(fill_halo) ) fill_halo_ = fill_halo
3842 if ( dim_type ==
'ZXY' &
3843 .OR. dim_type ==
'ZXHY' &
3844 .OR. dim_type ==
'ZXYH' &
3845 .OR. dim_type ==
'ZXHYH' )
then
3848 elseif ( dim_type ==
'ZHXY' &
3849 .OR. dim_type ==
'ZHXHY' &
3850 .OR. dim_type ==
'ZHXYH' )
then
3853 elseif( dim_type ==
'OXY' )
then
3856 elseif( dim_type ==
'LXY' )
then
3859 elseif( dim_type ==
'UXY' )
then
3863 log_error(
"FILE_CARTESC_write_var_3D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
3881 allocate( buf(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e) )
3883 if ( fill_halo_ )
then
3888 do k = dim1_s, dim1_e
3889 buf(k,i,j) = var(k,i,j)
3895 do j = dim3_s, dim3_e
3897 do k = dim1_s, dim1_e
3903 do j = dim3_s, dim3_e
3905 do k = dim1_s, dim1_e
3912 do i = dim2_s, dim2_e
3913 do k = dim1_s, dim1_e
3920 do i = dim2_s, dim2_e
3921 do k = dim1_s, dim1_e
3929 do j = dim3_s, dim3_e
3930 do i = dim2_s, dim2_e
3931 do k = dim1_s, dim1_e
3932 buf(k,i,j) = var(k,i,j)
3939 call file_write( vid, buf(:,:,:), nowdaysec, nowdaysec, start )
3977 integer,
intent(in) :: fid
3978 integer,
intent(in) :: vid
3979 real(RP),
intent(in) :: var(:,:,:)
3980 character(len=*),
intent(in) :: varname
3981 character(len=*),
intent(in) :: dim_type
3982 real(DP),
intent(in) :: timeintv
3984 integer,
intent(in),
optional :: timetarg
3985 real(DP),
intent(in),
optional :: timeofs
3986 logical,
intent(in),
optional :: fill_halo
3988 real(RP),
allocatable :: buf(:,:)
3990 integer :: dim1_S, dim1_E
3991 integer :: dim2_S, dim2_E
3993 real(DP) :: time_interval, nowtime
3997 logical :: fill_halo_
3998 real(DP) :: timeofs_
3999 integer :: rankidx(2)
4007 fill_halo_ = .false.
4008 if(
present(fill_halo) ) fill_halo_ = fill_halo
4011 if(
present(timeofs) ) timeofs_ = timeofs
4013 time_interval = timeintv
4014 step =
size(var(
isb,
jsb,:))
4019 if ( dim_type ==
'XYT' )
then
4032 log_error(
"FILE_CARTESC_write_var_3D_t",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
4040 allocate( buf(dim1_s:dim1_e,dim2_s:dim2_e) )
4042 if (
present(timetarg) )
then
4043 nowtime = timeofs_ + (timetarg-1) * time_interval
4047 if ( fill_halo_ )
then
4051 buf(i,j) = var(i,j,timetarg)
4056 do j = dim2_s, dim2_e
4062 do j = dim2_s, dim2_e
4069 do i = dim1_s, dim1_e
4075 do i = dim1_s, dim1_e
4082 do j = dim2_s, dim2_e
4083 do i = dim1_s, dim1_e
4084 buf(i,j) = var(i,j,timetarg)
4090 call file_write( vid, buf(:,:), nowtime, nowtime, start )
4099 if ( fill_halo_ )
then
4103 buf(i,j) = var(i,j,n)
4108 do j = dim2_s, dim2_e
4114 do j = dim2_s, dim2_e
4121 do i = dim1_s, dim1_e
4127 do i = dim1_s, dim1_e
4134 do j = dim2_s, dim2_e
4135 do i = dim1_s, dim1_e
4136 buf(i,j) = var(i,j,n)
4142 call file_write( vid, buf(:,:), nowtime, nowtime, start )
4145 nowtime = nowtime + time_interval
4184 integer,
intent(in) :: fid
4185 integer,
intent(in) :: vid
4186 real(RP),
intent(in) :: var(:,:,:,:)
4187 character(len=*),
intent(in) :: varname
4188 character(len=*),
intent(in) :: dim_type
4189 real(DP),
intent(in) :: timeintv
4191 integer,
intent(in),
optional :: timetarg
4192 real(DP),
intent(in),
optional :: timeofs
4193 logical,
intent(in),
optional :: fill_halo
4195 real(RP),
allocatable :: buf(:,:,:)
4197 integer :: dim1_S, dim1_E
4198 integer :: dim2_S, dim2_E
4199 integer :: dim3_S, dim3_E
4201 real(DP) :: time_interval, nowtime
4204 integer :: i, j, k, n
4205 logical :: fill_halo_
4206 real(DP) :: timeofs_
4207 integer :: rankidx(2)
4215 fill_halo_ = .false.
4216 if(
present(fill_halo) ) fill_halo_ = fill_halo
4219 if(
present(timeofs) ) timeofs_ = timeofs
4229 time_interval = timeintv
4244 if ( dim_type ==
'ZXYT' &
4245 .OR. dim_type ==
'ZXHYT' &
4246 .OR. dim_type ==
'ZXYHT' )
then
4249 elseif ( dim_type ==
'ZHXYT' )
then
4252 elseif( dim_type ==
'OXYT' )
then
4255 elseif( dim_type ==
'OHXYT' )
then
4258 elseif( dim_type ==
'LXYT' )
then
4261 elseif( dim_type ==
'LHXYT' )
then
4264 elseif( dim_type ==
'UXYT' )
then
4267 elseif( dim_type ==
'UHXYT' )
then
4271 log_error(
"FILE_CARTESC_write_var_4D",*)
'unsupported dimension type. Check! dim_type:', trim(dim_type),
', item:',trim(varname)
4275 allocate( buf(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e) )
4279 if (
present(timetarg) )
then
4280 nowtime = timeofs_ + (timetarg-1) * time_interval
4282 if ( fill_halo_ )
then
4286 do k = dim1_s, dim1_e
4287 buf(k,i,j) = var(k,i,j,timetarg)
4293 do j = dim3_s, dim3_e
4295 do k = dim1_s, dim1_e
4301 do j = dim3_s, dim3_e
4303 do k = dim1_s, dim1_e
4310 do i = dim2_s, dim2_e
4311 do k = dim1_s, dim1_e
4318 do i = dim2_s, dim2_e
4319 do k = dim1_s, dim1_e
4327 do j = dim3_s, dim3_e
4328 do i = dim2_s, dim2_e
4329 do k = dim1_s, dim1_e
4330 buf(k,i,j) = var(k,i,j,timetarg)
4337 call file_write( vid, buf(:,:,:), nowtime, nowtime, start )
4343 if ( fill_halo_ )
then
4347 do k = dim1_s, dim1_e
4348 buf(k,i,j) = var(k,i,j,n)
4354 do j = dim3_s, dim3_e
4356 do k = dim1_s, dim1_e
4362 do j = dim3_s, dim3_e
4364 do k = dim1_s, dim1_e
4371 do i = dim2_s, dim2_e
4372 do k = dim1_s, dim1_e
4379 do i = dim2_s, dim2_e
4380 do k = dim1_s, dim1_e
4388 do j = dim3_s, dim3_e
4389 do i = dim2_s, dim2_e
4390 do k = dim1_s, dim1_e
4391 buf(k,i,j) = var(k,i,j,n)
4398 call file_write( vid, buf(:,:,:), nowtime, nowtime, start )
4401 nowtime = nowtime + time_interval
4419 subroutine check_1d( &
4428 real(RP),
intent(in) :: expected(:)
4429 real(RP),
intent(in) :: buffer(:)
4430 character(len=*),
intent(in) :: name
4438 nmax =
size(expected)
4439 if (
size(buffer) /= nmax )
then
4440 log_error(
"check_1d",*)
'size of coordinate ('//trim(name)//
') is different:', nmax,
size(buffer)
4445 if ( abs(expected(n)) > eps )
then
4446 check = abs(buffer(n)-expected(n)) / abs(buffer(n)+expected(n)) * 2.0_rp
4448 check = abs(buffer(n)-expected(n))
4451 if ( check > file_cartesc_datacheck_criteria )
then
4452 log_error(
"check_1d",*)
'value of coordinate ('//trim(name)//
') at ', n,
' is different:', &
4453 expected(n), buffer(n), check
4459 end subroutine check_1d
4462 subroutine check_2d( &
4471 real(RP),
intent(in) :: expected(:,:)
4472 real(RP),
intent(in) :: buffer(:,:)
4473 character(len=*),
intent(in) :: name
4476 integer :: imax, jmax
4482 imax =
size(expected,1)
4483 jmax =
size(expected,2)
4484 if (
size(buffer,1) /= imax )
then
4485 log_error(
"check_2d",*)
'the first size of coordinate ('//trim(name)//
') is different:', imax,
size(buffer,1)
4488 if (
size(buffer,2) /= jmax )
then
4489 log_error(
"check_2d",*)
'the second size of coordinate ('//trim(name)//
') is different:', jmax,
size(buffer,2)
4495 if ( abs(expected(i,j)) > eps )
then
4496 check = abs(buffer(i,j)-expected(i,j)) / abs(buffer(i,j)+expected(i,j)) * 2.0_rp
4498 check = abs(buffer(i,j)-expected(i,j))
4501 if ( check > file_cartesc_datacheck_criteria )
then
4502 log_error(
"check_2d",*)
'value of coordinate ('//trim(name)//
') at (', i,
',', j,
') is different:', &
4503 expected(i,j), buffer(i,j), check
4510 end subroutine check_2d
4513 subroutine check_3d( &
4523 real(RP),
intent(in) :: expected(:,:,:)
4524 real(RP),
intent(in) :: buffer(:,:,:)
4525 character(len=*),
intent(in) :: name
4526 logical,
intent(in) :: transpose
4529 integer :: imax, jmax, kmax
4535 if ( transpose )
then
4536 kmax =
size(expected,3)
4537 imax =
size(expected,1)
4538 jmax =
size(expected,2)
4540 kmax =
size(expected,1)
4541 imax =
size(expected,2)
4542 jmax =
size(expected,3)
4544 if (
size(buffer,1) /= kmax )
then
4545 log_error(
"check_3d",*)
'the first size of coordinate ('//trim(name)//
') is different:', kmax,
size(buffer,1)
4548 if (
size(buffer,2) /= imax )
then
4549 log_error(
"check_3d",*)
'the second size of coordinate ('//trim(name)//
') is different:', imax,
size(buffer,2)
4552 if (
size(buffer,3) /= jmax )
then
4553 log_error(
"check_3d",*)
'the third size of coordinate ('//trim(name)//
') is different:', jmax,
size(buffer,3)
4557 if ( transpose )
then
4562 if ( abs(expected(k,i,j)) > eps )
then
4563 check = abs(buffer(i,j,k)-expected(k,i,j)) / abs(buffer(i,j,k)+expected(k,i,j)) * 2.0_rp
4565 check = abs(buffer(i,j,k)-expected(k,i,j))
4568 if ( check > file_cartesc_datacheck_criteria )
then
4569 log_error(
"check_3d",*)
'value of coordinate ('//trim(name)//
') at ', i,
',', j,
',', k,
' is different:', &
4570 expected(k,i,j), buffer(i,j,k), check
4580 if ( abs(expected(k,i,j)) > eps )
then
4581 check = abs(buffer(k,i,j)-expected(k,i,j)) / abs(buffer(k,i,j)+expected(k,i,j)) * 2.0_rp
4583 check = abs(buffer(k,i,j)-expected(k,i,j))
4586 if ( check > file_cartesc_datacheck_criteria )
then
4587 log_error(
"check_3d",*)
'value of coordinate ('//trim(name)//
') at ', k,
',', i,
',', j,
' is different:', &
4588 expected(k,i,j), buffer(k,i,j), check
4597 end subroutine check_3d
4620 call set_dimension(
'XY', 2, (/
'x' ,
'y' /),
ia*
ja, .true., area=
'cell_area', location=
'face' )
4621 call set_dimension(
'UY', 2, (/
'xh',
'y ' /),
ia*
ja, .true., area=
'cell_area_uy', location=
'edge1' )
4622 call set_dimension(
'XV', 2, (/
'x ',
'yh' /),
ia*
ja, .true., area=
'cell_area_xv', location=
'edge2' )
4630 area=
'cell_area', area_x=
'cell_area_zxy_x', area_y=
'cell_area_zxy_y', &
4631 volume=
'cell_volume', location=
'face' )
4633 area=
'cell_area', area_x=
'cell_area_wxy_x', area_y=
'cell_area_wxy_y', &
4634 volume=
'cell_volume_wxy', location=
'face' )
4636 area=
'cell_area_uy', area_x=
'cell_area_zuy_x', &
4637 volume=
'cell_volume_zuy', location=
'edge1' )
4639 area=
'cell_area_xv', area_y=
'cell_area_zxv_y', &
4640 volume=
'cell_volume_zxv', location=
'edge2' )
4642 area_x=
'cell_area_zuv_x', area_y=
'cell_area_zuv_y', &
4645 area_x=
'cell_area_wuy_x', &
4648 area_y=
'cell_area_wxv_y', &
4651 if (
okmax > 0 )
then
4652 call set_dimension(
'OXY', 3, (/
'oz',
'x ',
'y ' /),
okmax*
ia*
ja, .true., area=
'cell_area', volume=
'cell_volume_oxy', location=
'face', grid=
'ocean' )
4653 call set_dimension(
'OHXY', 3, (/
'ozh',
'x ',
'y ' /), (
okmax+1)*
ia*
ja, .true., area=
'cell_area', volume=
'cell_volume_oxy', location=
'face', grid=
'ocean' )
4655 if (
lkmax > 0 )
then
4656 call set_dimension(
'LXY', 3, (/
'lz',
'x ',
'y ' /),
lkmax*
ia*
ja, .true., area=
'cell_area', volume=
'cell_volume_lxy', location=
'face', grid=
'land' )
4657 call set_dimension(
'LHXY', 3, (/
'lzh',
'x ',
'y ' /), (
lkmax+1)*
ia*
ja, .true., area=
'cell_area', volume=
'cell_volume_lxy', location=
'face', grid=
'land' )
4659 if (
ukmax > 0 )
then
4660 call set_dimension(
'UXY', 3, (/
'uz',
'x ',
'y ' /),
ukmax*
ia*
ja, .true., area=
'cell_area', volume=
'cell_volume_uxy', location=
'face', grid=
'urban' )
4661 call set_dimension(
'UHXY', 3, (/
'uzh',
'x ',
'y ' /), (
ukmax+1)*
ia*
ja, .true., area=
'cell_area', volume=
'cell_volume_uxy', location=
'face', grid=
'urban' )
4668 file_cartesc_axis_info(1)%periodic = .true.
4669 file_cartesc_axis_info(2)%periodic = .true.
4671 file_cartesc_axis_info(1)%periodic = .false.
4672 file_cartesc_axis_info(2)%periodic = .false.
4676 file_cartesc_axis_info(3)%periodic = .true.
4677 file_cartesc_axis_info(4)%periodic = .true.
4679 file_cartesc_axis_info(3)%periodic = .false.
4680 file_cartesc_axis_info(4)%periodic = .false.
4686 file_cartesc_axis_info(1)%size_global (1) =
imax *
prc_num_x
4687 file_cartesc_axis_info(1)%start_global(1) =
is_ing -
ihalo
4688 file_cartesc_axis_info(1)%halo_global (1) = 0
4689 file_cartesc_axis_info(1)%halo_global (2) = 0
4690 file_cartesc_axis_info(1)%halo_local (1) = 0
4691 file_cartesc_axis_info(1)%halo_local (2) = 0
4693 file_cartesc_axis_info(1)%size_global (1) =
iag
4694 file_cartesc_axis_info(1)%start_global(1) =
isga
4695 file_cartesc_axis_info(1)%halo_global (1) =
ihalo
4696 file_cartesc_axis_info(1)%halo_global (2) =
ihalo
4697 file_cartesc_axis_info(1)%halo_local (1) =
ihalo
4698 file_cartesc_axis_info(1)%halo_local (2) =
ihalo
4699 if(
prc_has_w ) file_cartesc_axis_info(1)%halo_local(1) = 0
4700 if(
prc_has_e ) file_cartesc_axis_info(1)%halo_local(2) = 0
4703 file_cartesc_axis_info(2) = file_cartesc_axis_info(1)
4707 file_cartesc_axis_info(3)%size_global (1) =
jmax *
prc_num_y
4708 file_cartesc_axis_info(3)%start_global(1) =
js_ing -
jhalo
4709 file_cartesc_axis_info(3)%halo_global (1) = 0
4710 file_cartesc_axis_info(3)%halo_global (2) = 0
4711 file_cartesc_axis_info(3)%halo_local (1) = 0
4712 file_cartesc_axis_info(3)%halo_local (2) = 0
4714 file_cartesc_axis_info(3)%size_global (1) =
jag
4715 file_cartesc_axis_info(3)%start_global(1) =
jsga
4716 file_cartesc_axis_info(3)%halo_global (1) =
jhalo
4717 file_cartesc_axis_info(3)%halo_global (2) =
jhalo
4718 file_cartesc_axis_info(3)%halo_local (1) =
jhalo
4719 file_cartesc_axis_info(3)%halo_local (2) =
jhalo
4720 if(
prc_has_s ) file_cartesc_axis_info(3)%halo_local(1) = 0
4721 if(
prc_has_n ) file_cartesc_axis_info(3)%halo_local(2) = 0
4724 file_cartesc_axis_info(4) = file_cartesc_axis_info(3)
4730 subroutine set_dimension( name, ndims, dims, size, mapping, area, area_x, area_y, volume, location, grid )
4733 character(len=*),
intent(in) :: name
4734 integer,
intent(in) :: ndims
4735 character(len=*),
intent(in) :: dims(ndims)
4736 integer,
intent(in) :: size
4737 logical,
intent(in),
optional :: mapping
4738 character(len=*),
intent(in),
optional :: area
4739 character(len=*),
intent(in),
optional :: area_x
4740 character(len=*),
intent(in),
optional :: area_y
4741 character(len=*),
intent(in),
optional :: volume
4742 character(len=*),
intent(in),
optional :: location
4743 character(len=*),
intent(in),
optional :: grid
4745 integer,
save :: dimid = 0
4751 if ( dimid > file_cartesc_ndims )
then
4752 log_error(
"set_dimension",*)
'number of dimensions exceeds the limit', dimid, file_cartesc_ndims
4757 file_cartesc_dims(dimid)%name = name
4759 file_cartesc_dims(dimid)%name = trim(name)//
"T"
4761 file_cartesc_dims(dimid)%ndims = ndims
4762 file_cartesc_dims(dimid)%dims(1:ndims) = dims(:)
4763 file_cartesc_dims(dimid)%size =
size
4765 if (
present(mapping) )
then
4766 file_cartesc_dims(dimid)%mapping = mapping
4768 file_cartesc_dims(dimid)%mapping = .false.
4771 if (
present(area) )
then
4772 file_cartesc_dims(dimid)%area = area
4774 file_cartesc_dims(dimid)%area =
''
4776 if (
present(area_x) )
then
4777 file_cartesc_dims(dimid)%area = area_x
4779 file_cartesc_dims(dimid)%area =
''
4781 if (
present(area_y) )
then
4782 file_cartesc_dims(dimid)%area = area_y
4784 file_cartesc_dims(dimid)%area =
''
4787 if (
present(volume) )
then
4788 file_cartesc_dims(dimid)%volume = volume
4790 file_cartesc_dims(dimid)%volume =
''
4793 if (
present(location) )
then
4794 file_cartesc_dims(dimid)%location = location
4795 if (
present(grid) )
then
4796 file_cartesc_dims(dimid)%grid =
'grid_'//trim(grid)
4798 file_cartesc_dims(dimid)%grid =
'grid'
4801 file_cartesc_dims(dimid)%location =
''
4802 file_cartesc_dims(dimid)%grid =
''
4819 integer :: err, order
4820 integer :: sizes(3), subsizes(3), sub_off(3)
4823 order = mpi_order_fortran
4825 centertypexy = mpi_datatype_null
4826 centertypezx = mpi_datatype_null
4827 centertypezxy = mpi_datatype_null
4828 centertypezhxy = mpi_datatype_null
4829 centertypeocean = mpi_datatype_null
4830 centertypeland = mpi_datatype_null
4831 centertypeurban = mpi_datatype_null
4835 if(
rp == 8 ) etype = mpi_double
4844 startzxy(2:3) = startxy(1:2)
4846 countzxy(2:3) = countxy(1:2)
4857 call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypezxy, err)
4858 call mpi_type_commit(centertypezxy, err)
4862 startzhxy(2:3) = startxy(1:2)
4863 countzhxy(1) =
kmax+1
4864 countzhxy(2:3) = countxy(1:2)
4869 subsizes(1) =
kmax+1
4875 call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypezhxy, err)
4876 call mpi_type_commit(centertypezhxy, err)
4878 if (
okmax > 0 )
then
4881 startocean(2:3) = startxy(1:2)
4882 countocean(1) =
okmax
4883 countocean(2:3) = countxy(1:2)
4887 sub_off(1) =
oks - 1
4888 call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypeocean, err)
4889 call mpi_type_commit(centertypeocean, err)
4892 if (
lkmax > 0 )
then
4895 startland(2:3) = startxy(1:2)
4896 countland(1) =
lkmax
4897 countland(2:3) = countxy(1:2)
4901 sub_off(1) =
lks - 1
4902 call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypeland, err)
4903 call mpi_type_commit(centertypeland, err)
4906 if (
ukmax > 0 )
then
4909 starturban(2:3) = startxy(1:2)
4910 counturban(1) =
ukmax
4911 counturban(2:3) = countxy(1:2)
4915 sub_off(1) =
uks - 1
4916 call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypeurban, err)
4917 call mpi_type_commit(centertypeurban, err)
4921 startzx(1) =
khalo+1
4931 sub_off(2) =
isb - 1
4932 call mpi_type_create_subarray(2, sizes, subsizes, sub_off, order, etype, centertypezx, err)
4933 call mpi_type_commit(centertypezx, err)
4947 if( centertypexy /= mpi_datatype_null )
call mpi_type_free(centertypexy, err)
4948 if( centertypezx /= mpi_datatype_null )
call mpi_type_free(centertypezx, err)
4949 if( centertypezxy /= mpi_datatype_null )
call mpi_type_free(centertypezxy, err)
4950 if( centertypezhxy /= mpi_datatype_null )
call mpi_type_free(centertypezhxy, err)
4951 if( centertypeocean /= mpi_datatype_null )
call mpi_type_free(centertypeocean, err)
4952 if( centertypeland /= mpi_datatype_null )
call mpi_type_free(centertypeland, err)
4953 if( centertypeurban /= mpi_datatype_null )
call mpi_type_free(centertypeurban, err)