SCALE-RM
scale_fileio.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public procedure
31  !
32  public :: fileio_setup
33  public :: fileio_set_coordinates
34  public :: fileio_set_axes
35  public :: fileio_read
36  public :: fileio_write
37 
38  public :: fileio_create
39  public :: fileio_def_var
40  public :: fileio_enddef
41  public :: fileio_write_var
42  public :: fileio_close
43 
44  interface fileio_read
45  module procedure fileio_read_1d
46  module procedure fileio_read_2d
47  module procedure fileio_read_3d
48  module procedure fileio_read_4d
49  end interface fileio_read
50 
51  interface fileio_write
52  module procedure fileio_write_1d
53  module procedure fileio_write_2d
54  module procedure fileio_write_3d
55  module procedure fileio_write_3d_t
56  module procedure fileio_write_4d
57  end interface fileio_write
58 
59  interface fileio_write_var
60  module procedure fileio_write_var_1d
61  module procedure fileio_write_var_2d
62  module procedure fileio_write_var_3d
63  module procedure fileio_write_var_3d_t
64  module procedure fileio_write_var_4d
65  end interface fileio_write_var
66 
67  !-----------------------------------------------------------------------------
68  !-----------------------------------------------------------------------------
69  !
70  !++ Public parameters & variables
71  !
72  !-----------------------------------------------------------------------------
73  !
74  !++ Private procedure
75  !
76  !-----------------------------------------------------------------------------
77  !
78  !++ Private parameters & variables
79  !
80  real(RP), private, allocatable :: axis_lon (:,:) ! [deg]
81  real(RP), private, allocatable :: axis_lonx(:,:) ! [deg]
82  real(RP), private, allocatable :: axis_lony(:,:) ! [deg]
83  real(RP), private, allocatable :: axis_lonxy(:,:) ! [deg]
84  real(RP), private, allocatable :: axis_lat (:,:) ! [deg]
85  real(RP), private, allocatable :: axis_latx(:,:) ! [deg]
86  real(RP), private, allocatable :: axis_laty(:,:) ! [deg]
87  real(RP), private, allocatable :: axis_latxy(:,:) ! [deg]
88 
89  integer, private, parameter :: file_nfile_max = 512 ! number limit of file
90  ! Keep consistency with "FILE_MAX" in gtool_netcdf.c
91  logical, private, save :: file_axes_written(0:file_nfile_max-1)
92  ! whether axes have been written
93  ! ! fid starts from zero so index should start from zero
94  logical, private, save :: file_nozcoord(0:file_nfile_max-1)
95  ! whether nozcoord is true or false
96 
97  !-----------------------------------------------------------------------------
98 contains
99  !-----------------------------------------------------------------------------
101  subroutine fileio_setup
102  implicit none
103  !---------------------------------------------------------------------------
104 
105  if( io_l ) write(io_fid_log,*)
106  if( io_l ) write(io_fid_log,*) '++++++ Module[FIELIO] / Categ[ATMOS-RM IO] / Origin[SCALElib]'
107  if( io_l ) write(io_fid_log,*) '*** No namelists.'
108  if( io_l ) write(io_fid_log,*)
109  if( io_l ) write(io_fid_log,*) '*** NetCDF header information ***'
110  if( io_l ) write(io_fid_log,*) '*** Data source : ', trim(h_source)
111  if( io_l ) write(io_fid_log,*) '*** Institute : ', trim(h_institute)
112 
113  allocate( axis_lon(imaxb,jmaxb) )
114  allocate( axis_lonx(imaxb,jmaxb) )
115  allocate( axis_lony(imaxb,jmaxb) )
116  allocate( axis_lonxy(imaxb,jmaxb) )
117  allocate( axis_lat(imaxb,jmaxb) )
118  allocate( axis_latx(imaxb,jmaxb) )
119  allocate( axis_laty(imaxb,jmaxb) )
120  allocate( axis_latxy(imaxb,jmaxb) )
121 
122  return
123  end subroutine fileio_setup
124 
125  !-----------------------------------------------------------------------------
127  subroutine fileio_set_coordinates( &
128  LON, &
129  LONX, &
130  LONY, &
131  LONXY, &
132  LAT, &
133  LATX, &
134  LATY, &
135  LATXY, &
136  CZ, &
137  FZ )
138  use scale_const, only: &
139  d2r => const_d2r
140  implicit none
141 
142  real(RP), intent(in) :: LON (ia,ja)
143  real(RP), intent(in) :: LONX(ia,ja)
144  real(RP), intent(in) :: LONY(ia,ja)
145  real(RP), intent(in) :: LONXY(ia,ja)
146  real(RP), intent(in) :: LAT (ia,ja)
147  real(RP), intent(in) :: LATX(ia,ja)
148  real(RP), intent(in) :: LATY(ia,ja)
149  real(RP), intent(in) :: LATXY(ia,ja)
150  real(RP), intent(in) :: CZ ( ka,ia,ja)
151  real(RP), intent(in) :: FZ (0:ka,ia,ja)
152  !---------------------------------------------------------------------------
153 
154  axis_lon(:,:) = lon(isb:ieb,jsb:jeb) / d2r
155  axis_lonx(:,:) = lonx(isb:ieb,jsb:jeb) / d2r
156  axis_lony(:,:) = lony(isb:ieb,jsb:jeb) / d2r
157  axis_lonxy(:,:) = lonxy(isb:ieb,jsb:jeb) / d2r
158  axis_lat(:,:) = lat(isb:ieb,jsb:jeb) / d2r
159  axis_latx(:,:) = latx(isb:ieb,jsb:jeb) / d2r
160  axis_laty(:,:) = laty(isb:ieb,jsb:jeb) / d2r
161  axis_latxy(:,:) = latxy(isb:ieb,jsb:jeb) / d2r
162 
163  return
164  end subroutine fileio_set_coordinates
165 
166  !-----------------------------------------------------------------------------
168  subroutine fileio_set_axes( &
169  fid, &
170  dtype, &
171  xy )
172  use gtool_file, only: &
173  fileputaxis, &
174  filesettattr, &
175  fileputassociatedcoordinates
176  use scale_grid, only: &
177  grid_cz, &
178  grid_cx, &
179  grid_cy, &
180  grid_fz, &
181  grid_fx, &
182  grid_fy, &
183  grid_cdz, &
184  grid_cdx, &
185  grid_cdy, &
186  grid_fdz, &
187  grid_fdx, &
188  grid_fdy, &
189  grid_cbfz, &
190  grid_cbfx, &
191  grid_cbfy, &
192  grid_fbfz, &
193  grid_fbfx, &
194  grid_fbfy, &
195  grid_cxg, &
196  grid_cyg, &
197  grid_fxg, &
198  grid_fyg, &
199  grid_cbfxg, &
200  grid_cbfyg, &
201  grid_fbfxg, &
202  grid_fbfyg
203  use scale_land_grid, only: &
204  grid_lcz, &
205  grid_lfz
206  use scale_urban_grid, only: &
207  grid_ucz, &
208  grid_ufz
209  implicit none
210 
211  integer, intent(in) :: fid
212  integer, intent(in) :: dtype
213  logical, intent(in), optional :: xy
214 
215  character(len=2) :: AXIS_name(2)
216  logical :: xy_
217  !---------------------------------------------------------------------------
218 
219  if ( present(xy) ) then
220  xy_ = xy
221  else
222  xy_ = .false.
223  end if
224 
225  if ( .NOT. xy_ ) then
226  call fileputaxis( fid, 'z', 'Z', 'm', 'z', dtype, grid_cz(ks:ke) )
227  end if
228  call fileputaxis( fid, 'x', 'X', 'm', 'x', dtype, grid_cx(isb:ieb) )
229  call fileputaxis( fid, 'y', 'Y', 'm', 'y', dtype, grid_cy(jsb:jeb) )
230  if ( .NOT. xy_ ) then
231  call fileputaxis( fid, 'zh', 'Z (half level)', 'm', 'zh', dtype, grid_fz(ks:ke) )
232  end if
233  call fileputaxis( fid, 'xh', 'X (half level)', 'm', 'xh', dtype, grid_fx(isb:ieb) )
234  call fileputaxis( fid, 'yh', 'Y (half level)', 'm', 'yh', dtype, grid_fy(jsb:jeb) )
235 
236  if ( .NOT. xy_ ) then
237  call fileputaxis( fid, 'lz', 'LZ', 'm', 'lz', dtype, grid_lcz(lks:lke) )
238  call fileputaxis( fid, 'lzh', 'LZ (half level)', 'm', 'lzh', dtype, grid_lfz(lks:lke) )
239  call fileputaxis( fid, 'uz', 'UZ', 'm', 'uz', dtype, grid_ucz(uks:uke) )
240  call fileputaxis( fid, 'uzh', 'UZ (half level)', 'm', 'uzh', dtype, grid_ufz(uks:uke) )
241  end if
242 
243  if ( .NOT. xy_ ) then
244  call fileputaxis( fid, 'CZ', 'Atmos Grid Center Position Z', 'm', 'CZ', dtype, grid_cz )
245  end if
246  call fileputaxis( fid, 'CX', 'Atmos Grid Center Position X', 'm', 'CX', dtype, grid_cx )
247  call fileputaxis( fid, 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', dtype, grid_cy )
248  if ( .NOT. xy_ ) then
249  call fileputaxis( fid, 'FZ', 'Atmos Grid Face Position Z', 'm', 'FZ', dtype, grid_fz )
250  end if
251  call fileputaxis( fid, 'FX', 'Atmos Grid Face Position X', 'm', 'FX', dtype, grid_fx )
252  call fileputaxis( fid, 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', dtype, grid_fy )
253 
254  if ( .NOT. xy_ ) then
255  call fileputaxis( fid, 'CDZ', 'Grid Cell length Z', 'm', 'CZ', dtype, grid_cdz )
256  end if
257  call fileputaxis( fid, 'CDX', 'Grid Cell length X', 'm', 'CX', dtype, grid_cdx )
258  call fileputaxis( fid, 'CDY', 'Grid Cell length Y', 'm', 'CY', dtype, grid_cdy )
259  if ( .NOT. xy_ ) then
260  call fileputaxis( fid, 'FDZ', 'Grid distance Z', 'm', 'FDZ', dtype, grid_fdz )
261  end if
262  call fileputaxis( fid, 'FDX', 'Grid distance X', 'm', 'FDX', dtype, grid_fdx )
263  call fileputaxis( fid, 'FDY', 'Grid distance Y', 'm', 'FDY', dtype, grid_fdy )
264 
265  if ( .NOT. xy_ ) then
266  call fileputaxis( fid, 'LCZ', 'Land Grid Center Position Z', 'm', 'LCZ', dtype, grid_lcz )
267  call fileputaxis( fid, 'LFZ', 'Land Grid Face Position Z', 'm', 'LFZ', dtype, grid_lfz )
268  call fileputaxis( fid, 'LCDZ', 'Land Grid Cell length Z', 'm', 'LCZ', dtype, grid_lcz )
269 
270  call fileputaxis( fid, 'UCZ', 'Urban Grid Center Position Z', 'm', 'UCZ', dtype, grid_ucz )
271  call fileputaxis( fid, 'UFZ', 'Urban Grid Face Position Z', 'm', 'UFZ', dtype, grid_ufz )
272  call fileputaxis( fid, 'UCDZ', 'Urban Grid Cell length Z', 'm', 'UCZ', dtype, grid_ucz )
273  end if
274 
275  if ( .NOT. xy_ ) then
276  call fileputaxis( fid, 'CBFZ', 'Boundary factor Center Z', '1', 'CZ', dtype, grid_cbfz )
277  end if
278  call fileputaxis( fid, 'CBFX', 'Boundary factor Center X', '1', 'CX', dtype, grid_cbfx )
279  call fileputaxis( fid, 'CBFY', 'Boundary factor Center Y', '1', 'CY', dtype, grid_cbfy )
280  if ( .NOT. xy_ ) then
281  call fileputaxis( fid, 'FBFZ', 'Boundary factor Face Z', '1', 'CZ', dtype, grid_fbfz )
282  end if
283  call fileputaxis( fid, 'FBFX', 'Boundary factor Face X', '1', 'CX', dtype, grid_fbfx )
284  call fileputaxis( fid, 'FBFY', 'Boundary factor Face Y', '1', 'CY', dtype, grid_fbfy )
285 
286  call fileputaxis( fid, 'CXG', 'Grid Center Position X (global)', 'm', 'CXG', dtype, grid_cxg )
287  call fileputaxis( fid, 'CYG', 'Grid Center Position Y (global)', 'm', 'CYG', dtype, grid_cyg )
288  call fileputaxis( fid, 'FXG', 'Grid Face Position X (global)', 'm', 'FXG', dtype, grid_fxg )
289  call fileputaxis( fid, 'FYG', 'Grid Face Position Y (global)', 'm', 'FYG', dtype, grid_fyg )
290 
291  call fileputaxis( fid, 'CBFXG', 'Boundary factor Center X (global)', '1', 'CXG', dtype, grid_cbfxg )
292  call fileputaxis( fid, 'CBFYG', 'Boundary factor Center Y (global)', '1', 'CYG', dtype, grid_cbfyg )
293  call fileputaxis( fid, 'FBFXG', 'Boundary factor Face X (global)', '1', 'CXG', dtype, grid_fbfxg )
294  call fileputaxis( fid, 'FBFYG', 'Boundary factor Face Y (global)', '1', 'CYG', dtype, grid_fbfyg )
295 
296  ! associate coordinates
297  axis_name = (/'x ','y '/)
298  call fileputassociatedcoordinates( fid, 'lon' , 'longitude' , &
299  'degrees_east' , axis_name, dtype, axis_lon(:,:) )
300  axis_name = (/'xh','y '/)
301  call fileputassociatedcoordinates( fid, 'lon_uy', 'longitude (half level uy)', &
302  'degrees_east' , axis_name, dtype, axis_lonx(:,:) )
303  axis_name = (/'x ','yh'/)
304  call fileputassociatedcoordinates( fid, 'lon_xv', 'longitude (half level xv)', &
305  'degrees_east' , axis_name, dtype, axis_lony(:,:) )
306  axis_name = (/'xh','yh'/)
307  call fileputassociatedcoordinates( fid, 'lon_uv', 'longitude (half level uv)', &
308  'degrees_east' , axis_name, dtype, axis_lonxy(:,:) )
309  axis_name = (/'x ','y '/)
310  call fileputassociatedcoordinates( fid, 'lat' , 'latitude' , &
311  'degrees_north', axis_name, dtype, axis_lat(:,:) )
312  axis_name = (/'xh','y '/)
313  call fileputassociatedcoordinates( fid, 'lat_uy', 'latitude (half level uy)' , &
314  'degrees_north', axis_name, dtype, axis_latx(:,:) )
315  axis_name = (/'x ','yh'/)
316  call fileputassociatedcoordinates( fid, 'lat_xv', 'latitude (half level xv)' , &
317  'degrees_north', axis_name, dtype, axis_laty(:,:) )
318  axis_name = (/'xh','yh'/)
319  call fileputassociatedcoordinates( fid, 'lat_uv', 'latitude (half level uv)' , &
320  'degrees_north', axis_name, dtype, axis_latxy(:,:) )
321 
322  ! attributes
323  if ( .NOT. xy_ ) then
324  call filesettattr( fid, 'lz', 'positive', 'down' )
325  call filesettattr( fid, 'lzh', 'positive', 'down' )
326  call filesettattr( fid, 'uz', 'positive', 'down' )
327  call filesettattr( fid, 'uzh', 'positive', 'down' )
328  call filesettattr( fid, 'LCZ', 'positive', 'down' )
329  call filesettattr( fid, 'LFZ', 'positive', 'down' )
330  call filesettattr( fid, 'UCZ', 'positive', 'down' )
331  call filesettattr( fid, 'UFZ', 'positive', 'down' )
332  end if
333 
334  return
335  end subroutine fileio_set_axes
336 
337  !-----------------------------------------------------------------------------
339  subroutine fileio_read_1d( &
340  var, &
341  basename, &
342  varname, &
343  axistype, &
344  step )
345  use gtool_file, only: &
346  fileread
347  use scale_process, only: &
348  prc_myrank, &
350  implicit none
351 
352  real(RP), intent(out) :: var(:)
353  character(len=*), intent(in) :: basename
354  character(len=*), intent(in) :: varname
355  character(len=*), intent(in) :: axistype
356  integer, intent(in) :: step
357 
358  integer :: dim1_max, dim1_S, dim1_E
359  real(RP), allocatable :: var1D(:)
360  !---------------------------------------------------------------------------
361 
362  call prof_rapstart('FILE_I_NetCDF', 2)
363 
364  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 1D var: ', trim(varname)
365 
366  if ( axistype == 'Z' ) then
367  dim1_max = kmax
368  dim1_s = ks
369  dim1_e = ke
370  elseif( axistype == 'X' ) then
371  dim1_max = imaxb
372  dim1_s = isb
373  dim1_e = ieb
374  elseif( axistype == 'Y' ) then
375  dim1_max = jmaxb
376  dim1_s = jsb
377  dim1_e = jeb
378  else
379  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
380  call prc_mpistop
381  endif
382 
383  allocate( var1d(dim1_max) )
384 
385  call fileread( var1d(:), basename, varname, step, prc_myrank )
386  var(dim1_s:dim1_e) = var1d(1:dim1_max)
387 
388  deallocate( var1d )
389 
390  call prof_rapend ('FILE_I_NetCDF', 2)
391 
392  return
393  end subroutine fileio_read_1d
394 
395  !-----------------------------------------------------------------------------
397  subroutine fileio_read_2d( &
398  var, &
399  basename, &
400  varname, &
401  axistype, &
402  step )
403  use gtool_file, only: &
404  fileread
405  use scale_process, only: &
406  prc_myrank, &
408  implicit none
409 
410  real(RP), intent(out) :: var(:,:)
411  character(len=*), intent(in) :: basename
412  character(len=*), intent(in) :: varname
413  character(len=*), intent(in) :: axistype
414  integer, intent(in) :: step
415 
416  integer :: dim1_max, dim1_S, dim1_E
417  integer :: dim2_max, dim2_S, dim2_E
418  real(RP), allocatable :: var2D(:,:)
419  !---------------------------------------------------------------------------
420 
421  call prof_rapstart('FILE_I_NetCDF', 2)
422 
423  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 2D var: ', trim(varname)
424 
425  if ( axistype == 'XY' ) then
426  dim1_max = imaxb
427  dim2_max = jmaxb
428  dim1_s = isb
429  dim1_e = ieb
430  dim2_s = jsb
431  dim2_e = jeb
432  elseif( axistype == 'ZX' ) then
433  dim1_max = kmax
434  dim2_max = imaxb
435  dim1_s = ks
436  dim1_e = ke
437  dim2_s = isb
438  dim2_e = ieb
439  else
440  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
441  call prc_mpistop
442  endif
443 
444  allocate( var2d(dim1_max,dim2_max) )
445 
446  call fileread( var2d(:,:), basename, varname, step, prc_myrank )
447  var(dim1_s:dim1_e,dim2_s:dim2_e) = var2d(1:dim1_max,1:dim2_max)
448 
449  deallocate( var2d )
450 
451  call prof_rapend ('FILE_I_NetCDF', 2)
452 
453  return
454  end subroutine fileio_read_2d
455 
456  !-----------------------------------------------------------------------------
458  subroutine fileio_read_3d( &
459  var, &
460  basename, &
461  varname, &
462  axistype, &
463  step )
464  use gtool_file, only: &
465  fileread
466  use scale_process, only: &
467  prc_myrank, &
469  implicit none
470 
471  real(RP), intent(out) :: var(:,:,:)
472  character(len=*), intent(in) :: basename
473  character(len=*), intent(in) :: varname
474  character(len=*), intent(in) :: axistype
475  integer, intent(in) :: step
476 
477  integer :: dim1_max, dim1_S, dim1_E
478  integer :: dim2_max, dim2_S, dim2_E
479  integer :: dim3_max, dim3_S, dim3_E
480  real(RP), allocatable :: var3D(:,:,:)
481  !---------------------------------------------------------------------------
482 
483  call prof_rapstart('FILE_I_NetCDF', 2)
484 
485  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 3D var: ', trim(varname)
486 
487  if ( axistype == 'ZXY' ) then
488  dim1_max = kmax
489  dim2_max = imaxb
490  dim3_max = jmaxb
491  dim1_s = ks
492  dim1_e = ke
493  dim2_s = isb
494  dim2_e = ieb
495  dim3_s = jsb
496  dim3_e = jeb
497  else if ( axistype == 'XYT' ) then
498  dim1_max = imaxb
499  dim2_max = jmaxb
500  dim3_max = step
501  dim1_s = isb
502  dim1_e = ieb
503  dim2_s = jsb
504  dim2_e = jeb
505  dim3_s = 1
506  dim3_e = step
507  else if ( axistype == 'Land' ) then
508  dim1_max = lkmax
509  dim2_max = imaxb
510  dim3_max = jmaxb
511  dim1_s = lks
512  dim1_e = lke
513  dim2_s = isb
514  dim2_e = ieb
515  dim3_s = jsb
516  dim3_e = jeb
517  else if ( axistype == 'Urban' ) then
518  dim1_max = ukmax
519  dim2_max = imaxb
520  dim3_max = jmaxb
521  dim1_s = uks
522  dim1_e = uke
523  dim2_s = isb
524  dim2_e = ieb
525  dim3_s = jsb
526  dim3_e = jeb
527  else
528  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
529  call prc_mpistop
530  endif
531 
532  allocate( var3d(dim1_max,dim2_max,dim3_max) )
533 
534  call fileread( var3d(:,:,:), basename, varname, step, prc_myrank )
535  var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e) = var3d(1:dim1_max,1:dim2_max,1:dim3_max)
536 
537  deallocate( var3d )
538 
539  call prof_rapend ('FILE_I_NetCDF', 2)
540 
541  return
542  end subroutine fileio_read_3d
543 
544  !-----------------------------------------------------------------------------
546  subroutine fileio_read_4d( &
547  var, &
548  basename, &
549  varname, &
550  axistype, &
551  step )
552  use gtool_file, only: &
553  fileread
554  use scale_process, only: &
555  prc_myrank, &
557  implicit none
558 
559  real(RP), intent(out) :: var(:,:,:,:)
560  character(len=*), intent(in) :: basename
561  character(len=*), intent(in) :: varname
562  character(len=*), intent(in) :: axistype
563  integer, intent(in) :: step
564 
565  integer :: dim1_max, dim1_S, dim1_E
566  integer :: dim2_max, dim2_S, dim2_E
567  integer :: dim3_max, dim3_S, dim3_E
568  integer :: dim4_max, dim4_S, dim4_E
569  real(RP), allocatable :: var4D(:,:,:,:)
570  !---------------------------------------------------------------------------
571 
572  call prof_rapstart('FILE_I_NetCDF', 2)
573 
574  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 4D var: ', trim(varname)
575 
576  if ( axistype == 'ZXYT' ) then
577  dim1_max = kmax
578  dim2_max = imaxb
579  dim3_max = jmaxb
580  dim4_max = step
581  dim1_s = ks
582  dim1_e = ke
583  dim2_s = isb
584  dim2_e = ieb
585  dim3_s = jsb
586  dim3_e = jeb
587  dim4_s = 1
588  dim4_e = step
589  else
590  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
591  call prc_mpistop
592  endif
593 
594  allocate( var4d(dim1_max,dim2_max,dim3_max,dim4_max) )
595 
596  call fileread( var4d(:,:,:,:), basename, varname, step, prc_myrank )
597  var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e,dim4_s:dim4_e) = var4d(1:dim1_max,1:dim2_max,1:dim3_max,1:dim4_max)
598 
599  deallocate( var4d )
600 
601  call prof_rapend ('FILE_I_NetCDF', 2)
602 
603  return
604  end subroutine fileio_read_4d
605 
606  !-----------------------------------------------------------------------------
608  subroutine fileio_write_1d( &
609  var, &
610  basename, &
611  title, &
612  varname, &
613  desc, &
614  unit, &
615  axistype, &
616  datatype, &
617  date, &
618  subsec, &
619  append )
620  use gtool_file_h, only: &
621  file_real8, &
622  file_real4
623  use gtool_file, only: &
624  filecreate, &
625  fileaddvariable, &
626  filesetglobalattribute, &
627  filewrite
628  use scale_process, only: &
629  prc_masterrank, &
630  prc_myrank, &
632  use scale_rm_process, only: &
633  prc_2drank
634  use scale_time, only: &
635  nowdate => time_nowdate, &
636  nowms => time_nowms, &
637  nowsec => time_nowdaysec
638  implicit none
639 
640  real(RP), intent(in) :: var(:)
641  character(len=*), intent(in) :: basename
642  character(len=*), intent(in) :: title
643  character(len=*), intent(in) :: varname
644  character(len=*), intent(in) :: desc
645  character(len=*), intent(in) :: unit
646  character(len=*), intent(in) :: axistype
647  character(len=*), intent(in) :: datatype
648 
649  integer, optional, intent(in) :: date(6)
650  real(DP),optional, intent(in) :: subsec
651  logical, optional, intent(in) :: append
652 
653  integer :: dtype
654  character(len=2) :: dims(1)
655  integer :: dim1_max, dim1_S, dim1_E
656  real(RP), allocatable :: var1D(:)
657 
658  integer :: rankidx(2)
659  logical :: fileexisted
660  integer :: fid, vid
661  character(len=34) :: tunits
662  !---------------------------------------------------------------------------
663 
664  call prof_rapstart('FILE_O_NetCDF', 2)
665 
666  rankidx(1) = prc_2drank(prc_myrank,1)
667  rankidx(2) = prc_2drank(prc_myrank,2)
668 
669  if ( datatype == 'REAL8' ) then
670  dtype = file_real8
671  elseif( datatype == 'REAL4' ) then
672  dtype = file_real4
673  else
674  if ( rp == 8 ) then
675  dtype = file_real8
676  elseif( rp == 4 ) then
677  dtype = file_real4
678  else
679  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
680  call prc_mpistop
681  endif
682  endif
683 
684  call filecreate( fid, & ! [OUT]
685  fileexisted, & ! [OUT]
686  basename, & ! [IN]
687  title, & ! [IN]
688  h_source, & ! [IN]
689  h_institute, & ! [IN]
690  prc_masterrank, & ! [IN]
691  prc_myrank, & ! [IN]
692  rankidx, & ! [IN]
693  append = append ) ! [IN]
694 
695  if ( .NOT. fileexisted ) then ! only once
696  call fileio_set_axes( fid, dtype ) ! [IN]
697  if ( present( subsec ) ) then
698  call filesetglobalattribute( fid, "time", (/subsec/) )
699  else
700  call filesetglobalattribute( fid, "time", (/nowms/) )
701  end if
702  if ( present( date ) ) then
703  call getcftunits(tunits, date)
704  else
705  call getcftunits(tunits, nowdate)
706  end if
707  call filesetglobalattribute( fid, "time_units", tunits )
708  endif
709 
710  if ( axistype == 'Z' ) then
711  dims = (/'z'/)
712  dim1_max = kmax
713  dim1_s = ks
714  dim1_e = ke
715  elseif( axistype == 'X' ) then
716  dims = (/'x'/)
717  dim1_max = imaxb
718  dim1_s = isb
719  dim1_e = ieb
720  elseif( axistype == 'Y' ) then
721  dims = (/'y'/)
722  dim1_max = jmaxb
723  dim1_s = jsb
724  dim1_e = jeb
725  else
726  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
727  call prc_mpistop
728  endif
729 
730  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype ) ! [IN]
731 
732  allocate( var1d(dim1_max) )
733 
734  var1d(1:dim1_max) = var(dim1_s:dim1_e)
735  call filewrite( fid, vid, var1d(:), nowsec, nowsec ) ! [IN]
736 
737  deallocate( var1d )
738 
739  call prof_rapend ('FILE_O_NetCDF', 2)
740 
741  return
742  end subroutine fileio_write_1d
743 
744  !-----------------------------------------------------------------------------
746  subroutine fileio_write_2d( &
747  var, &
748  basename, &
749  title, &
750  varname, &
751  desc, &
752  unit, &
753  axistype, &
754  datatype, &
755  date, &
756  subsec, &
757  append, &
758  nohalo, &
759  nozcoord )
760  use gtool_file, only: &
761  rmiss
762  use gtool_file_h, only: &
763  file_real8, &
764  file_real4
765  use gtool_file, only: &
766  filecreate, &
767  fileaddvariable, &
768  filesetglobalattribute, &
769  filewrite
770  use scale_process, only: &
771  prc_masterrank, &
772  prc_myrank, &
774  use scale_rm_process, only: &
775  prc_2drank
776  use scale_time, only: &
777  nowdate => time_nowdate, &
778  nowms => time_nowms, &
779  nowsec => time_nowdaysec
780  implicit none
781 
782  real(RP), intent(in) :: var(:,:)
783  character(len=*), intent(in) :: basename
784  character(len=*), intent(in) :: title
785  character(len=*), intent(in) :: varname
786  character(len=*), intent(in) :: desc
787  character(len=*), intent(in) :: unit
788  character(len=*), intent(in) :: axistype
789  character(len=*), intent(in) :: datatype
790  integer, optional, intent(in) :: date(6)
791  real(DP),optional, intent(in) :: subsec
792  logical, optional, intent(in) :: append
793  logical, optional, intent(in) :: nohalo
794  logical, optional, intent(in) :: nozcoord
795 
796  real(RP) :: varhalo( size(var(:,1)), size(var(1,:)) )
797 
798  integer :: dtype
799  character(len=2) :: dims(2)
800  integer :: dim1_max, dim1_S, dim1_E
801  integer :: dim2_max, dim2_S, dim2_E
802  real(RP), allocatable :: var2D(:,:)
803 
804  integer :: rankidx(2)
805  logical :: fileexisted
806  integer :: fid, vid
807  integer :: i, j
808  logical :: nohalo_
809  character(len=34) :: tunits
810  !---------------------------------------------------------------------------
811 
812  call prof_rapstart('FILE_O_NetCDF', 2)
813 
814  nohalo_ = .false.
815  if ( present(nohalo) ) nohalo_ = nohalo
816 
817  rankidx(1) = prc_2drank(prc_myrank,1)
818  rankidx(2) = prc_2drank(prc_myrank,2)
819 
820  if ( datatype == 'REAL8' ) then
821  dtype = file_real8
822  elseif( datatype == 'REAL4' ) then
823  dtype = file_real4
824  else
825  if ( rp == 8 ) then
826  dtype = file_real8
827  elseif( rp == 4 ) then
828  dtype = file_real4
829  else
830  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
831  call prc_mpistop
832  endif
833  endif
834 
835  call filecreate( fid, & ! [OUT]
836  fileexisted, & ! [OUT]
837  basename, & ! [IN]
838  title, & ! [IN]
839  h_source, & ! [IN]
840  h_institute, & ! [IN]
841  prc_masterrank, & ! [IN]
842  prc_myrank, & ! [IN]
843  rankidx, & ! [IN]
844  append = append ) ! [IN]
845 
846  if ( axistype == 'XY' ) then
847  dims = (/'x','y'/)
848  dim1_max = imaxb
849  dim2_max = jmaxb
850  dim1_s = isb
851  dim1_e = ieb
852  dim2_s = jsb
853  dim2_e = jeb
854  elseif ( axistype == 'UY' ) then
855  dims = (/'xh','y '/)
856  dim1_max = imaxb
857  dim2_max = jmaxb
858  dim1_s = isb
859  dim1_e = ieb
860  dim2_s = jsb
861  dim2_e = jeb
862  elseif ( axistype == 'XV' ) then
863  dims = (/'x ','yh'/)
864  dim1_max = imaxb
865  dim2_max = jmaxb
866  dim1_s = isb
867  dim1_e = ieb
868  dim2_s = jsb
869  dim2_e = jeb
870  elseif ( axistype == 'UV' ) then
871  dims = (/'xh','yh'/)
872  dim1_max = imaxb
873  dim2_max = jmaxb
874  dim1_s = isb
875  dim1_e = ieb
876  dim2_s = jsb
877  dim2_e = jeb
878  elseif( axistype == 'ZX' ) then
879  dims = (/'z','x'/)
880  dim1_max = kmax
881  dim2_max = imaxb
882  dim1_s = ks
883  dim1_e = ke
884  dim2_s = isb
885  dim2_e = ieb
886  else
887  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
888  call prc_mpistop
889  endif
890 
891  if ( .NOT. fileexisted ) then ! only once
892  call fileio_set_axes( fid, dtype, nozcoord ) ! [IN]
893  if ( present( subsec ) ) then
894  call filesetglobalattribute( fid, "time", (/subsec/) )
895  else
896  call filesetglobalattribute( fid, "time", (/nowms/) )
897  end if
898  if ( present( date ) ) then
899  call getcftunits(tunits, date)
900  else
901  call getcftunits(tunits, nowdate)
902  end if
903  call filesetglobalattribute( fid, "time_units", tunits )
904  endif
905 
906  varhalo(:,:) = var(:,:)
907 
908  if ( nohalo_ ) then
909  ! W halo
910  do j = 1, ja
911  do i = 1, is-1
912  varhalo(i,j) = rmiss
913  end do
914  end do
915  ! E halo
916  do j = 1, ja
917  do i = ie+1, ia
918  varhalo(i,j) = rmiss
919  end do
920  end do
921  ! S halo
922  do j = 1, js-1
923  do i = 1, ia
924  varhalo(i,j) = rmiss
925  end do
926  end do
927  ! N halo
928  do j = je+1, ja
929  do i = 1, ia
930  varhalo(i,j) = rmiss
931  end do
932  end do
933  end if
934 
935  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype ) ! [IN]
936 
937  allocate( var2d(dim1_max,dim2_max) )
938 
939  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
940  call filewrite( fid, vid, var2d(:,:), nowsec, nowsec ) ! [IN]
941 
942  deallocate( var2d )
943 
944  call prof_rapend ('FILE_O_NetCDF', 2)
945 
946  return
947  end subroutine fileio_write_2d
948 
949  !-----------------------------------------------------------------------------
951  subroutine fileio_write_3d( &
952  var, &
953  basename, &
954  title, &
955  varname, &
956  desc, &
957  unit, &
958  axistype, &
959  datatype, &
960  date, &
961  subsec, &
962  append, &
963  nohalo )
964  use gtool_file, only: &
965  rmiss
966  use gtool_file_h, only: &
967  file_real8, &
968  file_real4
969  use gtool_file, only: &
970  filecreate, &
971  fileaddvariable, &
972  filesetglobalattribute, &
973  filewrite
974  use scale_process, only: &
975  prc_masterrank, &
976  prc_myrank, &
978  use scale_rm_process, only: &
979  prc_2drank
980  use scale_time, only: &
981  nowdate => time_nowdate, &
982  nowms => time_nowms, &
983  nowsec => time_nowdaysec
984  implicit none
985 
986  real(RP), intent(in) :: var(:,:,:)
987  character(len=*), intent(in) :: basename
988  character(len=*), intent(in) :: title
989  character(len=*), intent(in) :: varname
990  character(len=*), intent(in) :: desc
991  character(len=*), intent(in) :: unit
992  character(len=*), intent(in) :: axistype
993  character(len=*), intent(in) :: datatype
994  integer, optional, intent(in) :: date(6)
995  real(DP),optional, intent(in) :: subsec
996  logical, optional, intent(in) :: append
997  logical, optional, intent(in) :: nohalo
998 
999  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)), size(var(1,1,:)) )
1000 
1001  integer :: dtype
1002  character(len=2) :: dims(3)
1003  integer :: dim1_max, dim1_S, dim1_E
1004  integer :: dim2_max, dim2_S, dim2_E
1005  integer :: dim3_max, dim3_S, dim3_E
1006 
1007  real(RP), allocatable :: var3D(:,:,:)
1008 
1009  integer :: rankidx(2)
1010  logical :: append_sw
1011  logical :: fileexisted
1012  integer :: fid, vid
1013  integer :: i, j, k
1014  logical :: nohalo_
1015  character(len=34) :: tunits
1016  !---------------------------------------------------------------------------
1017 
1018  call prof_rapstart('FILE_O_NetCDF', 2)
1019 
1020  nohalo_ = .false.
1021  if ( present(nohalo) ) nohalo_ = nohalo
1022 
1023  rankidx(1) = prc_2drank(prc_myrank,1)
1024  rankidx(2) = prc_2drank(prc_myrank,2)
1025 
1026  if ( datatype == 'REAL8' ) then
1027  dtype = file_real8
1028  elseif( datatype == 'REAL4' ) then
1029  dtype = file_real4
1030  else
1031  if ( rp == 8 ) then
1032  dtype = file_real8
1033  elseif( rp == 4 ) then
1034  dtype = file_real4
1035  else
1036  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
1037  call prc_mpistop
1038  endif
1039  endif
1040 
1041  append_sw = .false.
1042  if ( present(append) ) then
1043  append_sw = append
1044  endif
1045 
1046  call filecreate( fid, & ! [OUT]
1047  fileexisted, & ! [OUT]
1048  basename, & ! [IN]
1049  title, & ! [IN]
1050  h_source, & ! [IN]
1051  h_institute, & ! [IN]
1052  prc_masterrank, & ! [IN]
1053  prc_myrank, & ! [IN]
1054  rankidx, & ! [IN]
1055  append = append_sw ) ! [IN]
1056 
1057  if ( .NOT. fileexisted ) then ! only once
1058  call fileio_set_axes( fid, dtype ) ! [IN]
1059  if ( present( subsec ) ) then
1060  call filesetglobalattribute( fid, "time", (/subsec/) )
1061  else
1062  call filesetglobalattribute( fid, "time", (/nowms/) )
1063  end if
1064  if ( present( date ) ) then
1065  call getcftunits(tunits, date)
1066  else
1067  call getcftunits(tunits, nowdate)
1068  end if
1069  call filesetglobalattribute( fid, "time_units", tunits )
1070  endif
1071 
1072  if ( axistype == 'ZXY' ) then
1073  dims = (/'z','x','y'/)
1074  dim1_max = kmax
1075  dim2_max = imaxb
1076  dim3_max = jmaxb
1077  dim1_s = ks
1078  dim1_e = ke
1079  dim2_s = isb
1080  dim2_e = ieb
1081  dim3_s = jsb
1082  dim3_e = jeb
1083  elseif( axistype == 'ZHXY' ) then
1084  dims = (/'zh','x ','y '/)
1085  dim1_max = kmax
1086  dim2_max = imaxb
1087  dim3_max = jmaxb
1088  dim1_s = ks
1089  dim1_e = ke
1090  dim2_s = isb
1091  dim2_e = ieb
1092  dim3_s = jsb
1093  dim3_e = jeb
1094  elseif( axistype == 'ZXHY' ) then
1095  dims = (/'z ','xh','y '/)
1096  dim1_max = kmax
1097  dim2_max = imaxb
1098  dim3_max = jmaxb
1099  dim1_s = ks
1100  dim1_e = ke
1101  dim2_s = isb
1102  dim2_e = ieb
1103  dim3_s = jsb
1104  dim3_e = jeb
1105  elseif( axistype == 'ZXYH' ) then
1106  dims = (/'z ','x ','yh'/)
1107  dim1_max = kmax
1108  dim2_max = imaxb
1109  dim3_max = jmaxb
1110  dim1_s = ks
1111  dim1_e = ke
1112  dim2_s = isb
1113  dim2_e = ieb
1114  dim3_s = jsb
1115  dim3_e = jeb
1116  elseif( axistype == 'Land' ) then
1117  dims = (/'lz','x ','y '/)
1118  dim1_max = lkmax
1119  dim2_max = imaxb
1120  dim3_max = jmaxb
1121  dim1_s = lks
1122  dim1_e = lke
1123  dim2_s = isb
1124  dim2_e = ieb
1125  dim3_s = jsb
1126  dim3_e = jeb
1127  elseif( axistype == 'Urban' ) then
1128  dims = (/'uz','x ','y '/)
1129  dim1_max = ukmax
1130  dim2_max = imaxb
1131  dim3_max = jmaxb
1132  dim1_s = uks
1133  dim1_e = uke
1134  dim2_s = isb
1135  dim2_e = ieb
1136  dim3_s = jsb
1137  dim3_e = jeb
1138  else
1139  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
1140  call prc_mpistop
1141  endif
1142 
1143  varhalo(:,:,:) = var(:,:,:)
1144 
1145  if ( nohalo_ ) then
1146  ! W halo
1147  do k = 1, dim1_max
1148  do j = 1, ja
1149  do i = 1, is-1
1150  varhalo(k,i,j) = rmiss
1151  end do
1152  end do
1153  end do
1154  ! E halo
1155  do k = 1, dim1_max
1156  do j = 1, ja
1157  do i = ie+1, ia
1158  varhalo(k,i,j) = rmiss
1159  end do
1160  end do
1161  end do
1162  ! S halo
1163  do k = 1, dim1_max
1164  do j = 1, js-1
1165  do i = 1, ia
1166  varhalo(k,i,j) = rmiss
1167  end do
1168  end do
1169  end do
1170  ! N halo
1171  do k = 1, dim1_max
1172  do j = je+1, ja
1173  do i = 1, ia
1174  varhalo(k,i,j) = rmiss
1175  end do
1176  end do
1177  end do
1178  end if
1179 
1180  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype ) ! [IN]
1181 
1182  allocate( var3d(dim1_max,dim2_max,dim3_max) )
1183 
1184  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
1185  call filewrite( fid, vid, var3d(:,:,:), nowsec, nowsec ) ! [IN]
1186 
1187  deallocate( var3d )
1188 
1189  call prof_rapend ('FILE_O_NetCDF', 2)
1190 
1191  return
1192  end subroutine fileio_write_3d
1193 
1194  !-----------------------------------------------------------------------------
1196  subroutine fileio_write_3d_t( &
1197  var, &
1198  basename, &
1199  title, &
1200  varname, &
1201  desc, &
1202  unit, &
1203  axistype, &
1204  datatype, &
1205  timeintv, &
1206  tsince, &
1207  append, &
1208  timetarg, &
1209  nohalo )
1210  use gtool_file, only: &
1211  rmiss
1212  use gtool_file_h, only: &
1213  file_real8, &
1214  file_real4
1215  use gtool_file, only: &
1216  filecreate, &
1217  fileputaxis, &
1218  fileaddvariable, &
1219  filewrite
1220  use scale_process, only: &
1221  prc_masterrank, &
1222  prc_myrank, &
1223  prc_mpistop
1224  use scale_rm_process, only: &
1225  prc_2drank
1226  implicit none
1227 
1228  real(RP), intent(in) :: var(:,:,:)
1229  character(len=*), intent(in) :: basename
1230  character(len=*), intent(in) :: title
1231  character(len=*), intent(in) :: varname
1232  character(len=*), intent(in) :: desc
1233  character(len=*), intent(in) :: unit
1234  character(len=*), intent(in) :: axistype
1235  character(len=*), intent(in) :: datatype
1236  real(RP), intent(in) :: timeintv
1237  integer , intent(in) :: tsince(6)
1238  logical, optional, intent(in) :: append
1239  integer, optional, intent(in) :: timetarg
1240  logical, optional, intent(in) :: nohalo
1241 
1242  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)) )
1243 
1244  integer :: dtype
1245  character(len=2) :: dims(2)
1246  integer :: dim1_max, dim1_S, dim1_E
1247  integer :: dim2_max, dim2_S, dim2_E
1248 
1249  real(RP), allocatable :: var2D(:,:)
1250  real(DP) :: time_interval, nowtime
1251 
1252  character(len=34) :: tunits
1253 
1254  integer :: rankidx(2)
1255  logical :: append_sw
1256  logical :: fileexisted
1257  integer :: fid, vid
1258  integer :: step
1259  integer :: i, j, n
1260  logical :: nohalo_
1261  !---------------------------------------------------------------------------
1262 
1263  call prof_rapstart('FILE_O_NetCDF', 2)
1264 
1265  nohalo_ = .false.
1266  if ( present(nohalo) ) nohalo_ = nohalo
1267 
1268  time_interval = timeintv
1269  step = size(var(isb,jsb,:))
1270 
1271  rankidx(1) = prc_2drank(prc_myrank,1)
1272  rankidx(2) = prc_2drank(prc_myrank,2)
1273 
1274  if ( datatype == 'REAL8' ) then
1275  dtype = file_real8
1276  elseif( datatype == 'REAL4' ) then
1277  dtype = file_real4
1278  else
1279  if ( rp == 8 ) then
1280  dtype = file_real8
1281  elseif( rp == 4 ) then
1282  dtype = file_real4
1283  else
1284  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
1285  call prc_mpistop
1286  endif
1287  endif
1288 
1289  append_sw = .false.
1290  if ( present(append) ) then
1291  append_sw = append
1292  endif
1293 
1294  write(tunits,'(a,i4.4,"-",i2.2,"-",i2.2," ",i2.2,":",i2.2,":",i2.2)') 'seconds since ', tsince
1295 
1296  call filecreate( fid, & ! [OUT]
1297  fileexisted, & ! [OUT]
1298  basename, & ! [IN]
1299  title, & ! [IN]
1300  h_source, & ! [IN]
1301  h_institute, & ! [IN]
1302  prc_masterrank, & ! [IN]
1303  prc_myrank, & ! [IN]
1304  rankidx, & ! [IN]
1305  time_units = tunits, & ! [IN]
1306  append = append_sw ) ! [IN]
1307 
1308  if ( .NOT. fileexisted ) then ! only once
1309  call fileio_set_axes( fid, dtype ) ! [IN]
1310  endif
1311 
1312  if ( axistype == 'XYT' ) then
1313  dims = (/'x','y'/)
1314  dim1_max = imaxb
1315  dim2_max = jmaxb
1316  dim1_s = isb
1317  dim1_e = ieb
1318  dim2_s = jsb
1319  dim2_e = jeb
1320  else
1321  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
1322  call prc_mpistop
1323  endif
1324 
1325  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype, time_interval ) ! [IN]
1326  allocate( var2d(dim1_max,dim2_max) )
1327 
1328  if ( present(timetarg) ) then
1329  varhalo(:,:) = var(:,:,timetarg)
1330 
1331  if ( nohalo_ ) then
1332  ! W halo
1333  do j = 1, ja
1334  do i = 1, is-1
1335  varhalo(i,j) = rmiss
1336  end do
1337  end do
1338  ! E halo
1339  do j = 1, ja
1340  do i = ie+1, ia
1341  varhalo(i,j) = rmiss
1342  end do
1343  end do
1344  ! S halo
1345  do j = 1, js-1
1346  do i = 1, ia
1347  varhalo(i,j) = rmiss
1348  end do
1349  end do
1350  ! N halo
1351  do j = je+1, ja
1352  do i = 1, ia
1353  varhalo(i,j) = rmiss
1354  end do
1355  end do
1356  end if
1357 
1358  nowtime = (timetarg-1) * time_interval
1359  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
1360  call filewrite( fid, vid, var2d(:,:), nowtime, nowtime ) ! [IN]
1361  else
1362  nowtime = 0.0_dp
1363  do n = 1, step
1364  varhalo(:,:) = var(:,:,n)
1365 
1366  if ( nohalo_ ) then
1367  ! W halo
1368  do j = 1, ja
1369  do i = 1, is-1
1370  varhalo(i,j) = rmiss
1371  end do
1372  end do
1373  ! E halo
1374  do j = 1, ja
1375  do i = ie+1, ia
1376  varhalo(i,j) = rmiss
1377  end do
1378  end do
1379  ! S halo
1380  do j = 1, js-1
1381  do i = 1, ia
1382  varhalo(i,j) = rmiss
1383  end do
1384  end do
1385  ! N halo
1386  do j = je+1, ja
1387  do i = 1, ia
1388  varhalo(i,j) = rmiss
1389  end do
1390  end do
1391  end if
1392 
1393  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
1394  call filewrite( fid, vid, var2d(:,:), nowtime, nowtime ) ! [IN]
1395  nowtime = nowtime + time_interval
1396  enddo
1397  endif
1398 
1399  deallocate( var2d )
1400 
1401  call prof_rapend ('FILE_O_NetCDF', 2)
1402 
1403  return
1404  end subroutine fileio_write_3d_t
1405 
1406  !-----------------------------------------------------------------------------
1408  subroutine fileio_write_4d( &
1409  var, &
1410  basename, &
1411  title, &
1412  varname, &
1413  desc, &
1414  unit, &
1415  axistype, &
1416  datatype, &
1417  timeintv, &
1418  tsince, &
1419  append, &
1420  timetarg, &
1421  nohalo )
1422  use gtool_file, only: &
1423  rmiss
1424  use gtool_file_h, only: &
1425  file_real8, &
1426  file_real4
1427  use gtool_file, only: &
1428  filecreate, &
1429  fileputaxis, &
1430  fileaddvariable, &
1431  filewrite
1432  use scale_process, only: &
1433  prc_masterrank, &
1434  prc_myrank, &
1435  prc_mpistop
1436  use scale_rm_process, only: &
1437  prc_2drank
1438  implicit none
1439 
1440  real(RP), intent(in) :: var(:,:,:,:)
1441  character(len=*), intent(in) :: basename
1442  character(len=*), intent(in) :: title
1443  character(len=*), intent(in) :: varname
1444  character(len=*), intent(in) :: desc
1445  character(len=*), intent(in) :: unit
1446  character(len=*), intent(in) :: axistype
1447  character(len=*), intent(in) :: datatype
1448  real(RP), intent(in) :: timeintv
1449  integer , intent(in) :: tsince(6)
1450  logical, optional, intent(in) :: append
1451  integer, optional, intent(in) :: timetarg
1452  logical, optional, intent(in) :: nohalo
1453 
1454  real(RP) :: varhalo( size(var(:,1,1,1)), size(var(1,:,1,1)), size(var(1,1,:,1)) )
1455 
1456  integer :: dtype
1457  character(len=2) :: dims(3)
1458  integer :: dim1_max, dim1_S, dim1_E
1459  integer :: dim2_max, dim2_S, dim2_E
1460  integer :: dim3_max, dim3_S, dim3_E
1461 
1462  real(RP), allocatable :: var3D(:,:,:)
1463  real(DP) :: time_interval, nowtime
1464 
1465  character(len=34) :: tunits
1466 
1467  integer :: rankidx(2)
1468  logical :: append_sw
1469  logical :: fileexisted
1470  integer :: fid, vid
1471  integer :: step
1472  integer :: i, j, k, n
1473  logical :: nohalo_
1474  !---------------------------------------------------------------------------
1475 
1476  call prof_rapstart('FILE_O_NetCDF', 2)
1477 
1478  nohalo_ = .false.
1479  if ( present(nohalo) ) nohalo_ = nohalo
1480 
1481  time_interval = timeintv
1482  step = size(var(ks,isb,jsb,:))
1483 
1484  rankidx(1) = prc_2drank(prc_myrank,1)
1485  rankidx(2) = prc_2drank(prc_myrank,2)
1486 
1487  if ( datatype == 'REAL8' ) then
1488  dtype = file_real8
1489  elseif( datatype == 'REAL4' ) then
1490  dtype = file_real4
1491  else
1492  if ( rp == 8 ) then
1493  dtype = file_real8
1494  elseif( rp == 4 ) then
1495  dtype = file_real4
1496  else
1497  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
1498  call prc_mpistop
1499  endif
1500  endif
1501 
1502  append_sw = .false.
1503  if ( present(append) ) then
1504  append_sw = append
1505  endif
1506 
1507  call getcftunits(tunits, tsince)
1508 
1509  call filecreate( fid, & ! [OUT]
1510  fileexisted, & ! [OUT]
1511  basename, & ! [IN]
1512  title, & ! [IN]
1513  h_source, & ! [IN]
1514  h_institute, & ! [IN]
1515  prc_masterrank, & ! [IN]
1516  prc_myrank, & ! [IN]
1517  rankidx, & ! [IN]
1518  time_units = tunits, & ! [IN]
1519  append = append_sw ) ! [IN]
1520 
1521  if ( .NOT. fileexisted ) then ! only once
1522  call fileio_set_axes( fid, dtype ) ! [IN]
1523  endif
1524 
1525  if ( axistype == 'ZXYT' ) then
1526  dims = (/'z','x','y'/)
1527  dim1_max = kmax
1528  dim2_max = imaxb
1529  dim3_max = jmaxb
1530  dim1_s = ks
1531  dim1_e = ke
1532  dim2_s = isb
1533  dim2_e = ieb
1534  dim3_s = jsb
1535  dim3_e = jeb
1536  else
1537  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
1538  call prc_mpistop
1539  endif
1540 
1541  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype, time_interval ) ! [IN]
1542  allocate( var3d(dim1_max,dim2_max,dim3_max) )
1543 
1544  if ( present(timetarg) ) then
1545  varhalo(:,:,:) = var(:,:,:,timetarg)
1546 
1547  if ( nohalo_ ) then
1548  ! W halo
1549  do k = 1, dim1_max
1550  do j = 1, ja
1551  do i = 1, is-1
1552  varhalo(k,i,j) = rmiss
1553  end do
1554  end do
1555  end do
1556  ! E halo
1557  do k = 1, dim1_max
1558  do j = 1, ja
1559  do i = ie+1, ia
1560  varhalo(k,i,j) = rmiss
1561  end do
1562  end do
1563  end do
1564  ! S halo
1565  do k = 1, dim1_max
1566  do j = 1, js-1
1567  do i = 1, ia
1568  varhalo(k,i,j) = rmiss
1569  end do
1570  end do
1571  end do
1572  ! N halo
1573  do k = 1, dim1_max
1574  do j = je+1, ja
1575  do i = 1, ia
1576  varhalo(k,i,j) = rmiss
1577  end do
1578  end do
1579  end do
1580  end if
1581 
1582  nowtime = (timetarg-1) * time_interval
1583  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
1584  call filewrite( fid, vid, var3d(:,:,:), nowtime, nowtime ) ! [IN]
1585  else
1586  nowtime = 0.0_dp
1587  do n = 1, step
1588  varhalo(:,:,:) = var(:,:,:,n)
1589 
1590  if ( nohalo_ ) then
1591  ! W halo
1592  do k = 1, dim1_max
1593  do j = 1, ja
1594  do i = 1, is-1
1595  varhalo(k,i,j) = rmiss
1596  end do
1597  end do
1598  end do
1599  ! E halo
1600  do k = 1, dim1_max
1601  do j = 1, ja
1602  do i = ie+1, ia
1603  varhalo(k,i,j) = rmiss
1604  end do
1605  end do
1606  end do
1607  ! S halo
1608  do k = 1, dim1_max
1609  do j = 1, js-1
1610  do i = 1, ia
1611  varhalo(k,i,j) = rmiss
1612  end do
1613  end do
1614  end do
1615  ! N halo
1616  do k = 1, dim1_max
1617  do j = je+1, ja
1618  do i = 1, ia
1619  varhalo(k,i,j) = rmiss
1620  end do
1621  end do
1622  end do
1623  end if
1624 
1625  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
1626  call filewrite( fid, vid, var3d(:,:,:), nowtime, nowtime ) ! [IN]
1627  nowtime = nowtime + time_interval
1628  enddo
1629  endif
1630 
1631  deallocate( var3d )
1632 
1633  call prof_rapend ('FILE_O_NetCDF', 2)
1634 
1635  return
1636  end subroutine fileio_write_4d
1637 
1638  ! private
1639 
1640  subroutine getcftunits(tunits, date)
1641  character(len=34), intent(out) :: tunits
1642  integer, intent(in) :: date(6)
1643 
1644  write(tunits,'(a,i4.4,"-",i2.2,"-",i2.2," ",i2.2,":",i2.2,":",i2.2)') 'seconds since ', date
1645 
1646  return
1647  end subroutine getcftunits
1648 
1649  !-----------------------------------------------------------------------------
1651  subroutine fileio_create( &
1652  fid, &
1653  basename, &
1654  title, &
1655  datatype, &
1656  date, &
1657  subsec, &
1658  append, &
1659  nozcoord )
1660  use gtool_file_h, only: &
1661  file_real8, &
1662  file_real4
1663  use gtool_file, only: &
1664  filecreate, &
1665  filesetglobalattribute
1666  use scale_process, only: &
1667  prc_masterrank, &
1668  prc_myrank, &
1669  prc_mpistop
1670  use scale_rm_process, only: &
1671  prc_2drank
1672  use scale_time, only: &
1673  nowdate => time_nowdate, &
1674  nowms => time_nowms, &
1675  nowsec => time_nowdaysec
1676  implicit none
1677 
1678  integer, intent(out) :: fid
1679  character(len=*), intent(in) :: basename
1680  character(len=*), intent(in) :: title
1681  character(len=*), intent(in) :: datatype
1682  integer, optional, intent(in) :: date(6)
1683  real(DP),optional, intent(in) :: subsec
1684  logical, optional, intent(in) :: append
1685  logical, optional, intent(in) :: nozcoord
1686 
1687  integer :: rankidx(2)
1688  integer :: dtype
1689  logical :: append_sw
1690  logical :: fileexisted
1691  character(len=34) :: tunits
1692  !---------------------------------------------------------------------------
1693 
1694  call prof_rapstart('FILE_O_NetCDF', 2)
1695 
1696  rankidx(1) = prc_2drank(prc_myrank,1)
1697  rankidx(2) = prc_2drank(prc_myrank,2)
1698 
1699  ! dtype is used to define the data type of axis variables in file
1700  if ( datatype == 'REAL8' ) then
1701  dtype = file_real8
1702  elseif( datatype == 'REAL4' ) then
1703  dtype = file_real4
1704  else
1705  if ( rp == 8 ) then
1706  dtype = file_real8
1707  elseif( rp == 4 ) then
1708  dtype = file_real4
1709  else
1710  write(*,*) 'xxx unsupported data type. Check!', trim(datatype)
1711  call prc_mpistop
1712  endif
1713  endif
1714 
1715  append_sw = .false.
1716  if ( present(append) ) then
1717  append_sw = append
1718  endif
1719 
1720  ! create a netCDF file if not already existed. Otherwise, open it.
1721  if ( present(date) ) then
1722  call getcftunits(tunits, date)
1723 
1724  call filecreate( fid, & ! [OUT]
1725  fileexisted, & ! [OUT]
1726  basename, & ! [IN]
1727  title, & ! [IN]
1728  h_source, & ! [IN]
1729  h_institute, & ! [IN]
1730  prc_masterrank, & ! [IN]
1731  prc_myrank, & ! [IN]
1732  rankidx, & ! [IN]
1733  time_units = tunits, & ! [IN]
1734  append = append_sw ) ! [IN]
1735  else
1736  call filecreate( fid, & ! [OUT]
1737  fileexisted, & ! [OUT]
1738  basename, & ! [IN]
1739  title, & ! [IN]
1740  h_source, & ! [IN]
1741  h_institute, & ! [IN]
1742  prc_masterrank, & ! [IN]
1743  prc_myrank, & ! [IN]
1744  rankidx, & ! [IN]
1745  append = append_sw ) ! [IN]
1746  endif
1747 
1748  if ( .NOT. fileexisted ) then ! do below only once when file is created
1749  call fileio_def_axes( fid, dtype, nozcoord ) ! [IN]
1750  file_axes_written(fid) = .false.
1751  if ( present( nozcoord ) ) then
1752  file_nozcoord(fid) = nozcoord
1753  else
1754  file_nozcoord(fid) = .false.
1755  endif
1756 
1757  if ( present( subsec ) ) then
1758  call filesetglobalattribute( fid, "time", (/subsec/) )
1759  else
1760  call filesetglobalattribute( fid, "time", (/nowms/) )
1761  end if
1762  if ( present( date ) ) then
1763  call getcftunits(tunits, date)
1764  else
1765  call getcftunits(tunits, nowdate)
1766  end if
1767  call filesetglobalattribute( fid, "time_units", tunits )
1768  endif
1769 
1770  call prof_rapend ('FILE_O_NetCDF', 2)
1771 
1772  return
1773  end subroutine fileio_create
1774 
1775  !-----------------------------------------------------------------------------
1777  subroutine fileio_enddef( &
1778  fid)
1779  use gtool_file, only: &
1780  fileenddef
1781  implicit none
1782 
1783  integer, intent(in) :: fid
1784 
1785  !---------------------------------------------------------------------------
1786 
1787  call prof_rapstart('FILE_O_NetCDF', 2)
1788 
1789  call fileenddef( fid ) ! [IN]
1790 
1791  ! At first write axis variables
1792  if ( .NOT. file_axes_written(fid) ) then
1793  call fileio_write_axes(fid, file_nozcoord(fid))
1794  file_axes_written(fid) = .true.
1795  endif
1796 
1797  call prof_rapend ('FILE_O_NetCDF', 2)
1798 
1799  return
1800  end subroutine fileio_enddef
1801 
1802  !-----------------------------------------------------------------------------
1804  subroutine fileio_close( &
1805  fid)
1806  use gtool_file, only: &
1807  fileclose
1808  implicit none
1809 
1810  integer, intent(in) :: fid
1811 
1812  !---------------------------------------------------------------------------
1813 
1814  call prof_rapstart('FILE_O_NetCDF', 2)
1815 
1816  call fileclose( fid ) ! [IN]
1817 
1818  call prof_rapend ('FILE_O_NetCDF', 2)
1819 
1820  return
1821  end subroutine fileio_close
1822 
1823  !-----------------------------------------------------------------------------
1825  subroutine fileio_def_axes( &
1826  fid, &
1827  dtype, &
1828  xy )
1829  use gtool_file, only: &
1830  filedefaxis, &
1831  filesettattr, &
1833  use scale_rm_process, only: &
1834  prc_num_x, &
1835  prc_num_y
1836  implicit none
1837 
1838  integer, intent(in) :: fid
1839  integer, intent(in) :: dtype
1840  logical, intent(in), optional :: xy
1841 
1842  character(len=2) :: AXIS_name(2)
1843  logical :: xy_
1844  integer :: IAG, JAG
1845  !---------------------------------------------------------------------------
1846 
1847  ! array size (global domain)
1848  iag = ihalo + imax*prc_num_x + ihalo
1849  jag = jhalo + jmax*prc_num_y + jhalo
1850 
1851  if ( present(xy) ) then
1852  xy_ = xy
1853  else
1854  xy_ = .false.
1855  end if
1856 
1857  if ( .NOT. xy_ ) then
1858  call filedefaxis( fid, 'z', 'Z', 'm', 'z', dtype, ke-ks+1 )
1859  end if
1860  call filedefaxis( fid, 'x', 'X', 'm', 'x', dtype, ieb-isb+1 )
1861  call filedefaxis( fid, 'y', 'Y', 'm', 'y', dtype, jeb-jsb+1 )
1862  if ( .NOT. xy_ ) then
1863  call filedefaxis( fid, 'zh', 'Z (half level)', 'm', 'zh', dtype, ke-ks+1 )
1864  end if
1865  call filedefaxis( fid, 'xh', 'X (half level)', 'm', 'xh', dtype, ieb-isb+1 )
1866  call filedefaxis( fid, 'yh', 'Y (half level)', 'm', 'yh', dtype, jeb-jsb+1 )
1867 
1868  if ( .NOT. xy_ ) then
1869  call filedefaxis( fid, 'lz', 'LZ', 'm', 'lz', dtype, lke-lks+1 )
1870  call filedefaxis( fid, 'lzh', 'LZ (half level)', 'm', 'lzh', dtype, lke-lks+1 )
1871  call filedefaxis( fid, 'uz', 'UZ', 'm', 'uz', dtype, uke-uks+1 )
1872  call filedefaxis( fid, 'uzh', 'UZ (half level)', 'm', 'uzh', dtype, uke-uks+1 )
1873  end if
1874 
1875 
1876  if ( .NOT. xy_ ) then
1877  call filedefaxis( fid, 'CZ', 'Atmos Grid Center Position Z', 'm', 'CZ', dtype, ka )
1878  end if
1879  call filedefaxis( fid, 'CX', 'Atmos Grid Center Position X', 'm', 'CX', dtype, ia )
1880  call filedefaxis( fid, 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', dtype, ja )
1881  if ( .NOT. xy_ ) then
1882  call filedefaxis( fid, 'FZ', 'Atmos Grid Face Position Z', 'm', 'FZ', dtype, ka+1 )
1883  end if
1884  call filedefaxis( fid, 'FX', 'Atmos Grid Face Position X', 'm', 'FX', dtype, ia+1 )
1885  call filedefaxis( fid, 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', dtype, ja+1 )
1886 
1887  if ( .NOT. xy_ ) then
1888  call filedefaxis( fid, 'CDZ', 'Grid Cell length Z', 'm', 'CZ', dtype, ka )
1889  end if
1890  call filedefaxis( fid, 'CDX', 'Grid Cell length X', 'm', 'CX', dtype, ia )
1891  call filedefaxis( fid, 'CDY', 'Grid Cell length Y', 'm', 'CY', dtype, ja )
1892  if ( .NOT. xy_ ) then
1893  call filedefaxis( fid, 'FDZ', 'Grid distance Z', 'm', 'FDZ', dtype, ka-1 )
1894  end if
1895  call filedefaxis( fid, 'FDX', 'Grid distance X', 'm', 'FDX', dtype, ia-1 )
1896  call filedefaxis( fid, 'FDY', 'Grid distance Y', 'm', 'FDY', dtype, ja-1 )
1897 
1898  if ( .NOT. xy_ ) then
1899  call filedefaxis( fid, 'LCZ', 'Land Grid Center Position Z', 'm', 'LCZ', dtype, lke-lks+1 )
1900  call filedefaxis( fid, 'LFZ', 'Land Grid Face Position Z', 'm', 'LFZ', dtype, lke-lks+2 )
1901  call filedefaxis( fid, 'LCDZ', 'Land Grid Cell length Z', 'm', 'LCZ', dtype, lke-lks+1 )
1902 
1903  call filedefaxis( fid, 'UCZ', 'Urban Grid Center Position Z', 'm', 'UCZ', dtype, uke-uks+1 )
1904  call filedefaxis( fid, 'UFZ', 'Urban Grid Face Position Z', 'm', 'UFZ', dtype, uke-uks+2 )
1905  call filedefaxis( fid, 'UCDZ', 'Urban Grid Cell length Z', 'm', 'UCZ', dtype, uke-uks+1 )
1906  end if
1907 
1908  if ( .NOT. xy_ ) then
1909  call filedefaxis( fid, 'CBFZ', 'Boundary factor Center Z', '1', 'CZ', dtype, ka )
1910  end if
1911  call filedefaxis( fid, 'CBFX', 'Boundary factor Center X', '1', 'CX', dtype, ia )
1912  call filedefaxis( fid, 'CBFY', 'Boundary factor Center Y', '1', 'CY', dtype, ja )
1913  if ( .NOT. xy_ ) then
1914  call filedefaxis( fid, 'FBFZ', 'Boundary factor Face Z', '1', 'CZ', dtype, ka )
1915  end if
1916  call filedefaxis( fid, 'FBFX', 'Boundary factor Face X', '1', 'CX', dtype, ia )
1917  call filedefaxis( fid, 'FBFY', 'Boundary factor Face Y', '1', 'CY', dtype, ja )
1918 
1919  call filedefaxis( fid, 'CXG', 'Grid Center Position X (global)', 'm', 'CXG', dtype, iag )
1920  call filedefaxis( fid, 'CYG', 'Grid Center Position Y (global)', 'm', 'CYG', dtype, jag )
1921  call filedefaxis( fid, 'FXG', 'Grid Face Position X (global)', 'm', 'FXG', dtype, iag+1 )
1922  call filedefaxis( fid, 'FYG', 'Grid Face Position Y (global)', 'm', 'FYG', dtype, jag+1 )
1923 
1924  call filedefaxis( fid, 'CBFXG', 'Boundary factor Center X (global)', '1', 'CXG', dtype, iag )
1925  call filedefaxis( fid, 'CBFYG', 'Boundary factor Center Y (global)', '1', 'CYG', dtype, jag )
1926  call filedefaxis( fid, 'FBFXG', 'Boundary factor Face X (global)', '1', 'CXG', dtype, iag )
1927  call filedefaxis( fid, 'FBFYG', 'Boundary factor Face Y (global)', '1', 'CYG', dtype, jag )
1928 
1929  ! associate coordinates
1930  axis_name = (/'x ','y '/)
1931  call filedefassociatedcoordinates( fid, 'lon' , 'longitude', &
1932  'degrees_east' , axis_name, dtype )
1933  axis_name = (/'xh','y '/)
1934  call filedefassociatedcoordinates( fid, 'lon_uy', 'longitude (half level uy)', &
1935  'degrees_east' , axis_name, dtype )
1936  axis_name = (/'x ','yh'/)
1937  call filedefassociatedcoordinates( fid, 'lon_xv', 'longitude (half level xv)', &
1938  'degrees_east' , axis_name, dtype )
1939  axis_name = (/'xh','yh'/)
1940  call filedefassociatedcoordinates( fid, 'lon_uv', 'longitude (half level uv)', &
1941  'degrees_east' , axis_name, dtype )
1942  axis_name = (/'x ','y '/)
1943  call filedefassociatedcoordinates( fid, 'lat' , 'latitude', &
1944  'degrees_north', axis_name, dtype )
1945  axis_name = (/'xh','y '/)
1946  call filedefassociatedcoordinates( fid, 'lat_uy', 'latitude (half level uy)', &
1947  'degrees_north', axis_name, dtype )
1948  axis_name = (/'x ','yh'/)
1949  call filedefassociatedcoordinates( fid, 'lat_xv', 'latitude (half level xv)', &
1950  'degrees_north', axis_name, dtype )
1951  axis_name = (/'xh','yh'/)
1952  call filedefassociatedcoordinates( fid, 'lat_uv', 'latitude (half level uv)', &
1953  'degrees_north', axis_name, dtype )
1954 
1955  ! attributes
1956  if ( .NOT. xy_ ) then
1957  call filesettattr( fid, 'lz', 'positive', 'down' )
1958  call filesettattr( fid, 'lzh', 'positive', 'down' )
1959  call filesettattr( fid, 'uz', 'positive', 'down' )
1960  call filesettattr( fid, 'uzh', 'positive', 'down' )
1961  call filesettattr( fid, 'LCZ', 'positive', 'down' )
1962  call filesettattr( fid, 'LFZ', 'positive', 'down' )
1963  call filesettattr( fid, 'UCZ', 'positive', 'down' )
1964  call filesettattr( fid, 'UFZ', 'positive', 'down' )
1965  end if
1966 
1967  return
1968  end subroutine fileio_def_axes
1969 
1970  !-----------------------------------------------------------------------------
1972  subroutine fileio_write_axes( &
1973  fid, &
1974  xy )
1975  use gtool_file, only: &
1976  filewriteaxis, &
1977  filewriteassociatedcoordinates
1978  use scale_grid, only: &
1979  grid_cz, &
1980  grid_cx, &
1981  grid_cy, &
1982  grid_fz, &
1983  grid_fx, &
1984  grid_fy, &
1985  grid_cdz, &
1986  grid_cdx, &
1987  grid_cdy, &
1988  grid_fdz, &
1989  grid_fdx, &
1990  grid_fdy, &
1991  grid_cbfz, &
1992  grid_cbfx, &
1993  grid_cbfy, &
1994  grid_fbfz, &
1995  grid_fbfx, &
1996  grid_fbfy, &
1997  grid_cxg, &
1998  grid_cyg, &
1999  grid_fxg, &
2000  grid_fyg, &
2001  grid_cbfxg, &
2002  grid_cbfyg, &
2003  grid_fbfxg, &
2004  grid_fbfyg
2005  use scale_land_grid, only: &
2006  grid_lcz, &
2007  grid_lfz
2008  use scale_urban_grid, only: &
2009  grid_ucz, &
2010  grid_ufz
2011  implicit none
2012 
2013  integer, intent(in) :: fid
2014  logical, intent(in), optional :: xy
2015 
2016  logical :: xy_
2017  !---------------------------------------------------------------------------
2018 
2019  if ( present(xy) ) then
2020  xy_ = xy
2021  else
2022  xy_ = .false.
2023  end if
2024 
2025  if ( .NOT. xy_ ) then
2026  call filewriteaxis( fid, 'z', grid_cz(ks:ke) )
2027  end if
2028  call filewriteaxis( fid, 'x', grid_cx(isb:ieb) )
2029  call filewriteaxis( fid, 'y', grid_cy(jsb:jeb) )
2030  if ( .NOT. xy_ ) then
2031  call filewriteaxis( fid, 'zh', grid_fz(ks:ke) )
2032  end if
2033  call filewriteaxis( fid, 'xh', grid_fx(isb:ieb) )
2034  call filewriteaxis( fid, 'yh', grid_fy(jsb:jeb) )
2035 
2036  if ( .NOT. xy_ ) then
2037  call filewriteaxis( fid, 'lz', grid_lcz(lks:lke) )
2038  call filewriteaxis( fid, 'lzh', grid_lfz(lks:lke) )
2039  call filewriteaxis( fid, 'uz', grid_ucz(uks:uke) )
2040  call filewriteaxis( fid, 'uzh', grid_ufz(uks:uke) )
2041  end if
2042 
2043  if ( .NOT. xy_ ) then
2044  call filewriteaxis( fid, 'CZ', grid_cz )
2045  end if
2046  call filewriteaxis( fid, 'CX', grid_cx )
2047  call filewriteaxis( fid, 'CY', grid_cy )
2048  if ( .NOT. xy_ ) then
2049  call filewriteaxis( fid, 'FZ', grid_fz )
2050  end if
2051  call filewriteaxis( fid, 'FX', grid_fx )
2052  call filewriteaxis( fid, 'FY', grid_fy )
2053 
2054  if ( .NOT. xy_ ) then
2055  call filewriteaxis( fid, 'CDZ', grid_cdz )
2056  end if
2057  call filewriteaxis( fid, 'CDX', grid_cdx )
2058  call filewriteaxis( fid, 'CDY', grid_cdy )
2059  if ( .NOT. xy_ ) then
2060  call filewriteaxis( fid, 'FDZ', grid_fdz )
2061  end if
2062  call filewriteaxis( fid, 'FDX', grid_fdx )
2063  call filewriteaxis( fid, 'FDY', grid_fdy )
2064 
2065  if ( .NOT. xy_ ) then
2066  call filewriteaxis( fid, 'LCZ', grid_lcz )
2067  call filewriteaxis( fid, 'LFZ', grid_lfz )
2068  call filewriteaxis( fid, 'LCDZ', grid_lcz )
2069 
2070  call filewriteaxis( fid, 'UCZ', grid_ucz )
2071  call filewriteaxis( fid, 'UFZ', grid_ufz )
2072  call filewriteaxis( fid, 'UCDZ', grid_ucz )
2073  end if
2074 
2075  if ( .NOT. xy_ ) then
2076  call filewriteaxis( fid, 'CBFZ', grid_cbfz )
2077  end if
2078  call filewriteaxis( fid, 'CBFX', grid_cbfx )
2079  call filewriteaxis( fid, 'CBFY', grid_cbfy )
2080  if ( .NOT. xy_ ) then
2081  call filewriteaxis( fid, 'FBFZ', grid_fbfz )
2082  end if
2083  call filewriteaxis( fid, 'FBFX', grid_fbfx )
2084  call filewriteaxis( fid, 'FBFY', grid_fbfy )
2085 
2086  call filewriteaxis( fid, 'CXG', grid_cxg )
2087  call filewriteaxis( fid, 'CYG', grid_cyg )
2088  call filewriteaxis( fid, 'FXG', grid_fxg )
2089  call filewriteaxis( fid, 'FYG', grid_fyg )
2090 
2091  call filewriteaxis( fid, 'CBFXG', grid_cbfxg )
2092  call filewriteaxis( fid, 'CBFYG', grid_cbfyg )
2093  call filewriteaxis( fid, 'FBFXG', grid_fbfxg )
2094  call filewriteaxis( fid, 'FBFYG', grid_fbfyg )
2095 
2096 
2097  ! associate coordinates
2098  call filewriteassociatedcoordinates( fid, 'lon' , axis_lon(:,:) )
2099  call filewriteassociatedcoordinates( fid, 'lon_uy', axis_lonx(:,:) )
2100  call filewriteassociatedcoordinates( fid, 'lon_xv', axis_lony(:,:) )
2101  call filewriteassociatedcoordinates( fid, 'lon_uv', axis_lonxy(:,:) )
2102  call filewriteassociatedcoordinates( fid, 'lat' , axis_lat(:,:) )
2103  call filewriteassociatedcoordinates( fid, 'lat_uy', axis_latx(:,:) )
2104  call filewriteassociatedcoordinates( fid, 'lat_xv', axis_laty(:,:) )
2105  call filewriteassociatedcoordinates( fid, 'lat_uv', axis_latxy(:,:) )
2106 
2107  return
2108  end subroutine fileio_write_axes
2109 
2110  !-----------------------------------------------------------------------------
2112  subroutine fileio_def_var( &
2113  fid, &
2114  vid, &
2115  varname, &
2116  desc, &
2117  unit, &
2118  axistype, &
2119  datatype, &
2120  timeintv )
2121  use gtool_file_h, only: &
2122  file_real8, &
2123  file_real4
2124  use gtool_file, only: &
2126  use scale_process, only: &
2127  prc_mpistop
2128  implicit none
2129 
2130  integer, intent(in) :: fid
2131  integer, intent(out) :: vid
2132  character(len=*), intent(in) :: varname
2133  character(len=*), intent(in) :: desc
2134  character(len=*), intent(in) :: unit
2135  character(len=*), intent(in) :: axistype
2136  character(len=*), intent(in) :: datatype
2137  real(RP), optional, intent(in) :: timeintv
2138 
2139  integer :: dtype, ndims
2140  character(len=2) :: dims(3)
2141  real(DP) :: time_interval
2142  !---------------------------------------------------------------------------
2143 
2144  call prof_rapstart('FILE_O_NetCDF', 2)
2145 
2146  if ( datatype == 'REAL8' ) then
2147  dtype = file_real8
2148  elseif( datatype == 'REAL4' ) then
2149  dtype = file_real4
2150  else
2151  if ( rp == 8 ) then
2152  dtype = file_real8
2153  elseif( rp == 4 ) then
2154  dtype = file_real4
2155  else
2156  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
2157  call prc_mpistop
2158  endif
2159  endif
2160 
2161  if ( axistype == 'Z' ) then ! 1D variable
2162  ndims = 1
2163  dims(1) = 'z'
2164  elseif( axistype == 'X' ) then
2165  ndims = 1
2166  dims(1) = 'x'
2167  elseif( axistype == 'Y' ) then
2168  ndims = 1
2169  dims(1) = 'y'
2170  elseif ( axistype == 'XY' ) then ! 2D variable
2171  ndims = 2
2172  dims(1) = 'x'
2173  dims(2) = 'y'
2174  elseif ( axistype == 'UY' ) then
2175  ndims = 2
2176  dims(1) = 'xh'
2177  dims(2) = 'y'
2178  elseif ( axistype == 'XV' ) then
2179  ndims = 2
2180  dims(1) = 'x'
2181  dims(2) = 'yh'
2182  elseif ( axistype == 'UV' ) then
2183  ndims = 2
2184  dims(1) = 'xh'
2185  dims(2) = 'yh'
2186  elseif( axistype == 'ZX' ) then
2187  ndims = 2
2188  dims(1) = 'z'
2189  dims(2) = 'x'
2190  elseif ( axistype == 'ZXY' ) then ! 3D variable
2191  ndims = 3
2192  dims = (/'z','x','y'/)
2193  elseif( axistype == 'ZHXY' ) then
2194  ndims = 3
2195  dims = (/'zh','x ','y '/)
2196  elseif( axistype == 'ZXHY' ) then
2197  ndims = 3
2198  dims = (/'z ','xh','y '/)
2199  elseif( axistype == 'ZXYH' ) then
2200  ndims = 3
2201  dims = (/'z ','x ','yh'/)
2202  elseif( axistype == 'Land' ) then
2203  ndims = 3
2204  dims = (/'lz','x ','y '/)
2205  elseif( axistype == 'Urban' ) then
2206  ndims = 3
2207  dims = (/'uz','x ','y '/)
2208  elseif ( axistype == 'XYT' ) then ! 3D variable with time dimension
2209  ndims = 2
2210  dims(1) = 'x'
2211  dims(2) = 'y'
2212  elseif ( axistype == 'ZXYT' ) then ! 4D variable
2213  ndims = 3
2214  dims = (/'z','x','y'/)
2215  else
2216  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
2217  call prc_mpistop
2218  endif
2219 
2220  if ( present(timeintv) ) then ! 3D/4D variable with time dimension
2221  time_interval = timeintv
2222  call filedefinevariable( fid, vid, varname, desc, unit, ndims, dims, dtype, & ! [IN]
2223  tint=time_interval ) ! [IN]
2224  else
2225  call filedefinevariable( fid, vid, varname, desc, unit, ndims, dims, dtype ) ! [IN]
2226  endif
2227 
2228  call prof_rapend ('FILE_O_NetCDF', 2)
2229 
2230  return
2231  end subroutine fileio_def_var
2232 
2233  !-----------------------------------------------------------------------------
2235  subroutine fileio_write_var_1d( &
2236  fid, &
2237  vid, &
2238  var, &
2239  varname, &
2240  axistype )
2241  use gtool_file, only: &
2242  filewritevar
2243  use scale_process, only: &
2244  prc_mpistop
2245  use scale_time, only: &
2246  nowsec => time_nowdaysec
2247  implicit none
2248 
2249  integer, intent(in) :: fid
2250  integer, intent(in) :: vid
2251  real(RP), intent(in) :: var(:)
2252  character(len=*), intent(in) :: varname
2253  character(len=*), intent(in) :: axistype
2254 
2255  integer :: dim1_max, dim1_S, dim1_E
2256  real(RP), allocatable :: var1D(:)
2257 
2258  !---------------------------------------------------------------------------
2259 
2260  call prof_rapstart('FILE_O_NetCDF', 2)
2261 
2262  if ( axistype == 'Z' ) then
2263  dim1_max = kmax
2264  dim1_s = ks
2265  dim1_e = ke
2266  elseif( axistype == 'X' ) then
2267  dim1_max = imaxb
2268  dim1_s = isb
2269  dim1_e = ieb
2270  elseif( axistype == 'Y' ) then
2271  dim1_max = jmaxb
2272  dim1_s = jsb
2273  dim1_e = jeb
2274  else
2275  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
2276  call prc_mpistop
2277  endif
2278 
2279  allocate( var1d(dim1_max) )
2280 
2281  var1d(1:dim1_max) = var(dim1_s:dim1_e)
2282  call filewritevar( vid, var1d(:), nowsec, nowsec ) ! [IN]
2283 
2284  deallocate( var1d )
2285 
2286  call prof_rapend ('FILE_O_NetCDF', 2)
2287 
2288  return
2289  end subroutine fileio_write_var_1d
2290 
2291  !-----------------------------------------------------------------------------
2293  subroutine fileio_write_var_2d( &
2294  fid, &
2295  vid, &
2296  var, &
2297  varname, &
2298  axistype, &
2299  nohalo )
2300  use gtool_file, only: &
2301  rmiss
2302  use gtool_file, only: &
2303  filewritevar
2304  use scale_process, only: &
2305  prc_mpistop
2306  use scale_time, only: &
2307  nowsec => time_nowdaysec
2308  implicit none
2309 
2310  integer, intent(in) :: fid
2311  integer, intent(in) :: vid
2312  real(RP), intent(in) :: var(:,:)
2313  character(len=*), intent(in) :: varname
2314  character(len=*), intent(in) :: axistype
2315  logical, optional, intent(in) :: nohalo
2316 
2317  real(RP) :: varhalo( size(var(:,1)), size(var(1,:)) )
2318 
2319  integer :: dim1_max, dim1_S, dim1_E
2320  integer :: dim2_max, dim2_S, dim2_E
2321  real(RP), allocatable :: var2D(:,:)
2322 
2323  integer :: i, j
2324  logical :: nohalo_
2325  !---------------------------------------------------------------------------
2326 
2327  call prof_rapstart('FILE_O_NetCDF', 2)
2328 
2329  nohalo_ = .false.
2330  if ( present(nohalo) ) nohalo_ = nohalo
2331 
2332  if ( axistype == 'XY' ) then
2333  dim1_max = imaxb
2334  dim2_max = jmaxb
2335  dim1_s = isb
2336  dim1_e = ieb
2337  dim2_s = jsb
2338  dim2_e = jeb
2339  elseif ( axistype == 'UY' ) then
2340  dim1_max = imaxb
2341  dim2_max = jmaxb
2342  dim1_s = isb
2343  dim1_e = ieb
2344  dim2_s = jsb
2345  dim2_e = jeb
2346  elseif ( axistype == 'XV' ) then
2347  dim1_max = imaxb
2348  dim2_max = jmaxb
2349  dim1_s = isb
2350  dim1_e = ieb
2351  dim2_s = jsb
2352  dim2_e = jeb
2353  elseif ( axistype == 'UV' ) then
2354  dim1_max = imaxb
2355  dim2_max = jmaxb
2356  dim1_s = isb
2357  dim1_e = ieb
2358  dim2_s = jsb
2359  dim2_e = jeb
2360  elseif( axistype == 'ZX' ) then
2361  dim1_max = kmax
2362  dim2_max = imaxb
2363  dim1_s = ks
2364  dim1_e = ke
2365  dim2_s = isb
2366  dim2_e = ieb
2367  else
2368  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
2369  call prc_mpistop
2370  endif
2371 
2372  varhalo(:,:) = var(:,:)
2373 
2374  if ( nohalo_ ) then
2375  ! W halo
2376  do j = 1, ja
2377  do i = 1, is-1
2378  varhalo(i,j) = rmiss
2379  end do
2380  end do
2381  ! E halo
2382  do j = 1, ja
2383  do i = ie+1, ia
2384  varhalo(i,j) = rmiss
2385  end do
2386  end do
2387  ! S halo
2388  do j = 1, js-1
2389  do i = 1, ia
2390  varhalo(i,j) = rmiss
2391  end do
2392  end do
2393  ! N halo
2394  do j = je+1, ja
2395  do i = 1, ia
2396  varhalo(i,j) = rmiss
2397  end do
2398  end do
2399  end if
2400 
2401  allocate( var2d(dim1_max,dim2_max) )
2402 
2403  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
2404  call filewritevar( vid, var2d(:,:), nowsec, nowsec ) ! [IN]
2405 
2406  deallocate( var2d )
2407 
2408  call prof_rapend ('FILE_O_NetCDF', 2)
2409 
2410  return
2411  end subroutine fileio_write_var_2d
2412 
2413  !-----------------------------------------------------------------------------
2415  subroutine fileio_write_var_3d( &
2416  fid, &
2417  vid, &
2418  var, &
2419  varname, &
2420  axistype, &
2421  nohalo )
2422  use gtool_file, only: &
2423  rmiss
2424  use gtool_file, only: &
2425  filewritevar
2426  use scale_process, only: &
2427  prc_mpistop
2428  use scale_time, only: &
2429  nowsec => time_nowdaysec
2430  implicit none
2431 
2432  integer, intent(in) :: fid
2433  integer, intent(in) :: vid
2434  real(RP), intent(in) :: var(:,:,:)
2435  character(len=*), intent(in) :: varname
2436  character(len=*), intent(in) :: axistype
2437  logical, optional, intent(in) :: nohalo
2438 
2439  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)), size(var(1,1,:)) )
2440 
2441  integer :: dim1_max, dim1_S, dim1_E
2442  integer :: dim2_max, dim2_S, dim2_E
2443  integer :: dim3_max, dim3_S, dim3_E
2444 
2445  real(RP), allocatable :: var3D(:,:,:)
2446 
2447  integer :: i, j, k
2448  logical :: nohalo_
2449  !---------------------------------------------------------------------------
2450 
2451  call prof_rapstart('FILE_O_NetCDF', 2)
2452 
2453  nohalo_ = .false.
2454  if ( present(nohalo) ) nohalo_ = nohalo
2455 
2456  if ( axistype == 'ZXY' ) then
2457  dim1_max = kmax
2458  dim2_max = imaxb
2459  dim3_max = jmaxb
2460  dim1_s = ks
2461  dim1_e = ke
2462  dim2_s = isb
2463  dim2_e = ieb
2464  dim3_s = jsb
2465  dim3_e = jeb
2466  elseif( axistype == 'ZHXY' ) then
2467  dim1_max = kmax
2468  dim2_max = imaxb
2469  dim3_max = jmaxb
2470  dim1_s = ks
2471  dim1_e = ke
2472  dim2_s = isb
2473  dim2_e = ieb
2474  dim3_s = jsb
2475  dim3_e = jeb
2476  elseif( axistype == 'ZXHY' ) then
2477  dim1_max = kmax
2478  dim2_max = imaxb
2479  dim3_max = jmaxb
2480  dim1_s = ks
2481  dim1_e = ke
2482  dim2_s = isb
2483  dim2_e = ieb
2484  dim3_s = jsb
2485  dim3_e = jeb
2486  elseif( axistype == 'ZXYH' ) then
2487  dim1_max = kmax
2488  dim2_max = imaxb
2489  dim3_max = jmaxb
2490  dim1_s = ks
2491  dim1_e = ke
2492  dim2_s = isb
2493  dim2_e = ieb
2494  dim3_s = jsb
2495  dim3_e = jeb
2496  elseif( axistype == 'Land' ) then
2497  dim1_max = lkmax
2498  dim2_max = imaxb
2499  dim3_max = jmaxb
2500  dim1_s = lks
2501  dim1_e = lke
2502  dim2_s = isb
2503  dim2_e = ieb
2504  dim3_s = jsb
2505  dim3_e = jeb
2506  elseif( axistype == 'Urban' ) then
2507  dim1_max = ukmax
2508  dim2_max = imaxb
2509  dim3_max = jmaxb
2510  dim1_s = uks
2511  dim1_e = uke
2512  dim2_s = isb
2513  dim2_e = ieb
2514  dim3_s = jsb
2515  dim3_e = jeb
2516  else
2517  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
2518  call prc_mpistop
2519  endif
2520 
2521  varhalo(:,:,:) = var(:,:,:)
2522 
2523  if ( nohalo_ ) then
2524  ! W halo
2525  do k = 1, dim1_max
2526  do j = 1, ja
2527  do i = 1, is-1
2528  varhalo(k,i,j) = rmiss
2529  end do
2530  end do
2531  end do
2532  ! E halo
2533  do k = 1, dim1_max
2534  do j = 1, ja
2535  do i = ie+1, ia
2536  varhalo(k,i,j) = rmiss
2537  end do
2538  end do
2539  end do
2540  ! S halo
2541  do k = 1, dim1_max
2542  do j = 1, js-1
2543  do i = 1, ia
2544  varhalo(k,i,j) = rmiss
2545  end do
2546  end do
2547  end do
2548  ! N halo
2549  do k = 1, dim1_max
2550  do j = je+1, ja
2551  do i = 1, ia
2552  varhalo(k,i,j) = rmiss
2553  end do
2554  end do
2555  end do
2556  end if
2557 
2558  allocate( var3d(dim1_max,dim2_max,dim3_max) )
2559 
2560  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
2561  call filewritevar( vid, var3d(:,:,:), nowsec, nowsec ) ! [IN]
2562 
2563  deallocate( var3d )
2564 
2565  call prof_rapend ('FILE_O_NetCDF', 2)
2566 
2567  return
2568  end subroutine fileio_write_var_3d
2569 
2570  !-----------------------------------------------------------------------------
2572  subroutine fileio_write_var_3d_t( &
2573  fid, &
2574  vid, &
2575  var, &
2576  varname, &
2577  axistype, &
2578  timeintv, &
2579  timetarg, &
2580  nohalo )
2581  use gtool_file, only: &
2582  rmiss
2583  use gtool_file, only: &
2584  filewritevar
2585  use scale_process, only: &
2586  prc_mpistop
2587  implicit none
2588 
2589  integer, intent(in) :: fid
2590  integer, intent(in) :: vid
2591  real(RP), intent(in) :: var(:,:,:)
2592  character(len=*), intent(in) :: varname
2593  character(len=*), intent(in) :: axistype
2594  real(RP), intent(in) :: timeintv
2595  integer, optional, intent(in) :: timetarg
2596  logical, optional, intent(in) :: nohalo
2597 
2598  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)) )
2599 
2600  integer :: dim1_max, dim1_S, dim1_E
2601  integer :: dim2_max, dim2_S, dim2_E
2602 
2603  real(RP), allocatable :: var2D(:,:)
2604  real(DP) :: time_interval, nowtime
2605 
2606  integer :: step
2607  integer :: i, j, n
2608  logical :: nohalo_
2609  !---------------------------------------------------------------------------
2610 
2611  call prof_rapstart('FILE_O_NetCDF', 2)
2612 
2613  nohalo_ = .false.
2614  if ( present(nohalo) ) nohalo_ = nohalo
2615 
2616  time_interval = timeintv
2617  step = size(var(isb,jsb,:))
2618 
2619  if ( axistype == 'XYT' ) then
2620  dim1_max = imaxb
2621  dim2_max = jmaxb
2622  dim1_s = isb
2623  dim1_e = ieb
2624  dim2_s = jsb
2625  dim2_e = jeb
2626  else
2627  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
2628  call prc_mpistop
2629  endif
2630 
2631  allocate( var2d(dim1_max,dim2_max) )
2632 
2633  if ( present(timetarg) ) then
2634  varhalo(:,:) = var(:,:,timetarg)
2635 
2636  if ( nohalo_ ) then
2637  ! W halo
2638  do j = 1, ja
2639  do i = 1, is-1
2640  varhalo(i,j) = rmiss
2641  end do
2642  end do
2643  ! E halo
2644  do j = 1, ja
2645  do i = ie+1, ia
2646  varhalo(i,j) = rmiss
2647  end do
2648  end do
2649  ! S halo
2650  do j = 1, js-1
2651  do i = 1, ia
2652  varhalo(i,j) = rmiss
2653  end do
2654  end do
2655  ! N halo
2656  do j = je+1, ja
2657  do i = 1, ia
2658  varhalo(i,j) = rmiss
2659  end do
2660  end do
2661  end if
2662 
2663  nowtime = (timetarg-1) * time_interval
2664  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
2665  call filewritevar( vid, var2d(:,:), nowtime, nowtime ) ! [IN]
2666  else
2667  nowtime = 0.0_dp
2668  do n = 1, step
2669  varhalo(:,:) = var(:,:,n)
2670 
2671  if ( nohalo_ ) then
2672  ! W halo
2673  do j = 1, ja
2674  do i = 1, is-1
2675  varhalo(i,j) = rmiss
2676  end do
2677  end do
2678  ! E halo
2679  do j = 1, ja
2680  do i = ie+1, ia
2681  varhalo(i,j) = rmiss
2682  end do
2683  end do
2684  ! S halo
2685  do j = 1, js-1
2686  do i = 1, ia
2687  varhalo(i,j) = rmiss
2688  end do
2689  end do
2690  ! N halo
2691  do j = je+1, ja
2692  do i = 1, ia
2693  varhalo(i,j) = rmiss
2694  end do
2695  end do
2696  end if
2697 
2698  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
2699  call filewritevar( vid, var2d(:,:), nowtime, nowtime ) ! [IN]
2700  nowtime = nowtime + time_interval
2701  enddo
2702  endif
2703 
2704  deallocate( var2d )
2705 
2706  call prof_rapend ('FILE_O_NetCDF', 2)
2707 
2708  return
2709  end subroutine fileio_write_var_3d_t
2710 
2711  !-----------------------------------------------------------------------------
2713  subroutine fileio_write_var_4d( &
2714  fid, &
2715  vid, &
2716  var, &
2717  varname, &
2718  axistype, &
2719  timeintv, &
2720  timetarg, &
2721  nohalo )
2722  use gtool_file, only: &
2723  rmiss
2724  use gtool_file, only: &
2725  filewritevar
2726  use scale_process, only: &
2727  prc_mpistop
2728  implicit none
2729 
2730  integer, intent(in) :: fid
2731  integer, intent(in) :: vid
2732  real(RP), intent(in) :: var(:,:,:,:)
2733  character(len=*), intent(in) :: varname
2734  character(len=*), intent(in) :: axistype
2735  real(RP), intent(in) :: timeintv
2736  integer, optional, intent(in) :: timetarg
2737  logical, optional, intent(in) :: nohalo
2738 
2739  real(RP) :: varhalo( size(var(:,1,1,1)), size(var(1,:,1,1)), size(var(1,1,:,1)) )
2740 
2741  integer :: dim1_max, dim1_S, dim1_E
2742  integer :: dim2_max, dim2_S, dim2_E
2743  integer :: dim3_max, dim3_S, dim3_E
2744 
2745  real(RP), allocatable :: var3D(:,:,:)
2746  real(DP) :: time_interval, nowtime
2747 
2748 
2749  integer :: step
2750  integer :: i, j, k, n
2751  logical :: nohalo_
2752  !---------------------------------------------------------------------------
2753 
2754  call prof_rapstart('FILE_O_NetCDF', 2)
2755 
2756  nohalo_ = .false.
2757  if ( present(nohalo) ) nohalo_ = nohalo
2758 
2759  time_interval = timeintv
2760  step = size(var(ks,isb,jsb,:))
2761 
2762  if ( axistype == 'ZXYT' ) then
2763  dim1_max = kmax
2764  dim2_max = imaxb
2765  dim3_max = jmaxb
2766  dim1_s = ks
2767  dim1_e = ke
2768  dim2_s = isb
2769  dim2_e = ieb
2770  dim3_s = jsb
2771  dim3_e = jeb
2772  else
2773  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
2774  call prc_mpistop
2775  endif
2776 
2777  allocate( var3d(dim1_max,dim2_max,dim3_max) )
2778 
2779  if ( present(timetarg) ) then
2780  varhalo(:,:,:) = var(:,:,:,timetarg)
2781 
2782  if ( nohalo_ ) then
2783  ! W halo
2784  do k = 1, dim1_max
2785  do j = 1, ja
2786  do i = 1, is-1
2787  varhalo(k,i,j) = rmiss
2788  end do
2789  end do
2790  end do
2791  ! E halo
2792  do k = 1, dim1_max
2793  do j = 1, ja
2794  do i = ie+1, ia
2795  varhalo(k,i,j) = rmiss
2796  end do
2797  end do
2798  end do
2799  ! S halo
2800  do k = 1, dim1_max
2801  do j = 1, js-1
2802  do i = 1, ia
2803  varhalo(k,i,j) = rmiss
2804  end do
2805  end do
2806  end do
2807  ! N halo
2808  do k = 1, dim1_max
2809  do j = je+1, ja
2810  do i = 1, ia
2811  varhalo(k,i,j) = rmiss
2812  end do
2813  end do
2814  end do
2815  end if
2816 
2817  nowtime = (timetarg-1) * time_interval
2818  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
2819  call filewritevar( vid, var3d(:,:,:), nowtime, nowtime ) ! [IN]
2820  else
2821  nowtime = 0.0_dp
2822  do n = 1, step
2823  varhalo(:,:,:) = var(:,:,:,n)
2824 
2825  if ( nohalo_ ) then
2826  ! W halo
2827  do k = 1, dim1_max
2828  do j = 1, ja
2829  do i = 1, is-1
2830  varhalo(k,i,j) = rmiss
2831  end do
2832  end do
2833  end do
2834  ! E halo
2835  do k = 1, dim1_max
2836  do j = 1, ja
2837  do i = ie+1, ia
2838  varhalo(k,i,j) = rmiss
2839  end do
2840  end do
2841  end do
2842  ! S halo
2843  do k = 1, dim1_max
2844  do j = 1, js-1
2845  do i = 1, ia
2846  varhalo(k,i,j) = rmiss
2847  end do
2848  end do
2849  end do
2850  ! N halo
2851  do k = 1, dim1_max
2852  do j = je+1, ja
2853  do i = 1, ia
2854  varhalo(k,i,j) = rmiss
2855  end do
2856  end do
2857  end do
2858  end if
2859 
2860  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
2861  call filewritevar( vid, var3d(:,:,:), nowtime, nowtime ) ! [IN]
2862  nowtime = nowtime + time_interval
2863  enddo
2864  endif
2865 
2866  deallocate( var3d )
2867 
2868  call prof_rapend ('FILE_O_NetCDF', 2)
2869 
2870  return
2871  end subroutine fileio_write_var_4d
2872 
2873 end module scale_fileio
integer, public imax
of computational cells: x
integer, public prc_num_x
x length of 2D processor topology
integer, public is
start point of inner domain: x, local
subroutine, public filecreate(fid, existed, basename, title, source, institution, master, myrank, rankidx, single, time_units, append)
Definition: gtool_file.f90:181
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public je
end point of inner domain: y, local
subroutine fileio_read_1d(var, basename, varname, axistype, step)
Read 1D data from file.
real(rp), dimension(:), allocatable, public grid_cyg
center coordinate [m]: y, global
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:), allocatable, public grid_cbfyg
center buffer factor [0-1]: y, global
integer, public jeb
real(dp), public time_nowms
subsecond part of current time [millisec]
Definition: scale_time.F90:66
real(dp), parameter, public rmiss
Definition: gtool_file.f90:141
real(rp), dimension(:), allocatable, public grid_cxg
center coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fdy
y-length of grid(j+1) to grid(j) [m]
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine fileio_def_axes(fid, dtype, xy)
define axis variables in the file
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
subroutine fileio_write_var_3d(fid, vid, var, varname, axistype, nohalo)
Write 3D data to file.
real(rp), dimension(:), allocatable, public grid_fbfy
face buffer factor [0-1]: y
subroutine, public fileio_setup
Setup.
subroutine, public filedefinevariable(fid, vid, varname, desc, units, ndims, dims, dtype, tint, tavg)
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
module STDIO
Definition: scale_stdio.F90:12
subroutine, public fileenddef(fid)
integer, public ke
end point of inner domain: z, local
subroutine fileio_write_3d_t(var, basename, title, varname, desc, unit, axistype, datatype, timeintv, tsince, append, timetarg, nohalo)
Write 3D data with time dimension to file.
integer, public imaxb
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
subroutine fileio_write_4d(var, basename, title, varname, desc, unit, axistype, datatype, timeintv, tsince, append, timetarg, nohalo)
Write 4D data to file.
module GRID (cartesian) for land
subroutine fileio_write_axes(fid, xy)
write axis to the file
subroutine fileio_write_var_4d(fid, vid, var, varname, axistype, timeintv, timetarg, nohalo)
Write 4D data to file.
subroutine fileio_read_4d(var, basename, varname, axistype, step)
Read 4D data from file.
integer, public jmaxb
module FILE I/O (netcdf)
real(rp), dimension(:), allocatable, public grid_fx
face coordinate [m]: x, local
subroutine fileio_write_var_1d(fid, vid, var, varname, axistype)
Write 1D data to file.
real(rp), dimension(:), allocatable, public grid_lfz
face coordinate [m]: z, local=global
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
module grid index
subroutine, public fileio_def_var(fid, vid, varname, desc, unit, axistype, datatype, timeintv)
Define a variable to file.
integer, public ia
of x whole cells (local, with HALO)
character(len=h_mid), public h_source
for file header
Definition: scale_stdio.F90:50
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]
subroutine, public fileclose(fid)
real(rp), dimension(:), allocatable, public grid_ufz
face coordinate [m]: z, local=global
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 filedefassociatedcoordinates(fid, name, desc, units, dim_names, dtype)
Definition: gtool_file.f90:869
subroutine, public filesettattr(fid, vname, key, val)
subroutine fileio_write_2d(var, basename, title, varname, desc, unit, axistype, datatype, date, subsec, append, nohalo, nozcoord)
Write 2D data to file.
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 grid_fbfz
face buffer factor [0-1]: z
subroutine fileio_write_3d(var, basename, title, varname, desc, unit, axistype, datatype, date, subsec, append, nohalo)
Write 3D data to file.
subroutine fileio_write_var_3d_t(fid, vid, var, varname, axistype, timeintv, timetarg, nohalo)
Write 3D data with time dimension to file.
subroutine, public fileio_create(fid, basename, title, datatype, date, subsec, append, nozcoord)
Create/open a netCDF file.
integer, public jhalo
of halo cells: y
integer, public js
start point of inner domain: y, local
real(rp), dimension(:), allocatable, public grid_cbfx
center buffer factor [0-1]: x
module TIME
Definition: scale_time.F90:15
module PROCESS
real(rp), dimension(:), allocatable, public grid_fbfyg
face buffer factor [0-1]: y, global
real(rp), dimension(:), allocatable, public grid_ucz
center coordinate [m]: z, local=global
module CONSTANT
Definition: scale_const.F90:14
integer, parameter, public prc_masterrank
master process in each communicator
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:), allocatable, public grid_cbfz
center buffer factor [0-1]: z
integer, public prc_myrank
process num in local communicator
real(rp), dimension(:), allocatable, public grid_cx
center coordinate [m]: x, local
module GRID (cartesian)
subroutine, public fileio_enddef(fid)
Exit netCDF file define mode.
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
module RM PROCESS
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
subroutine getcftunits(tunits, date)
character(len=h_mid), public h_institute
for file header
Definition: scale_stdio.F90:51
module PRECISION
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
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
integer, parameter, public file_real4
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, public isb
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
subroutine, public fileio_close(fid)
Close a netCDF file.
module land grid index
subroutine, public fileio_set_coordinates(LON, LONX, LONY, LONXY, LAT, LATX, LATY, LATXY, CZ, FZ)
set latlon and z
module FILE I/O HEADER
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public fileio_set_axes(fid, dtype, xy)
write axis to the file
integer, public jsb
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
integer, parameter, public rp
real(rp), dimension(:), allocatable, public grid_cdx
x-length of control volume [m]
integer, public jmax
of computational cells: y
subroutine fileio_write_var_2d(fid, vid, var, varname, axistype, nohalo)
Write 2D data to file.
real(rp), dimension(:), allocatable, public grid_fyg
face coordinate [m]: y, global
real(rp), dimension(:), allocatable, public grid_cy
center coordinate [m]: y, local
integer, parameter, public file_real8
subroutine fileio_read_2d(var, basename, varname, axistype, step)
Read 2D data from file.
real(rp), dimension(:), allocatable, public grid_cbfxg
center buffer factor [0-1]: x, global
subroutine fileio_read_3d(var, basename, varname, axistype, step)
Read 3D data from file.
integer, public ihalo
of halo cells: x
subroutine fileio_write_1d(var, basename, title, varname, desc, unit, axistype, datatype, date, subsec, append)
Write 1D data to file.
subroutine, public filedefaxis(fid, name, desc, units, dim_name, dtype, dim_size)
Definition: gtool_file.f90:564
module urban grid index
integer, public ja
of y whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local