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  interface fileio_read
39  module procedure fileio_read_1d
40  module procedure fileio_read_2d
41  module procedure fileio_read_3d
42  module procedure fileio_read_4d
43  end interface fileio_read
44 
45  interface fileio_write
46  module procedure fileio_write_1d
47  module procedure fileio_write_2d
48  module procedure fileio_write_3d
49  module procedure fileio_write_3d_t
50  module procedure fileio_write_4d
51  end interface fileio_write
52 
53  !-----------------------------------------------------------------------------
54  !
55  !++ Public parameters & variables
56  !
57  !-----------------------------------------------------------------------------
58  !
59  !++ Private procedure
60  !
61  !-----------------------------------------------------------------------------
62  !
63  !++ Private parameters & variables
64  !
65  real(RP), private, allocatable :: axis_lon (:,:) ! [deg]
66  real(RP), private, allocatable :: axis_lonx(:,:) ! [deg]
67  real(RP), private, allocatable :: axis_lony(:,:) ! [deg]
68  real(RP), private, allocatable :: axis_lonxy(:,:) ! [deg]
69  real(RP), private, allocatable :: axis_lat (:,:) ! [deg]
70  real(RP), private, allocatable :: axis_latx(:,:) ! [deg]
71  real(RP), private, allocatable :: axis_laty(:,:) ! [deg]
72  real(RP), private, allocatable :: axis_latxy(:,:) ! [deg]
73  !-----------------------------------------------------------------------------
74 contains
75  !-----------------------------------------------------------------------------
77  subroutine fileio_setup
78  implicit none
79  !---------------------------------------------------------------------------
80 
81  if( io_l ) write(io_fid_log,*)
82  if( io_l ) write(io_fid_log,*) '++++++ Module[FIELIO] / Categ[ATMOS-RM IO] / Origin[SCALElib]'
83  if( io_l ) write(io_fid_log,*) '*** No namelists.'
84  if( io_l ) write(io_fid_log,*)
85  if( io_l ) write(io_fid_log,*) '*** NetCDF header information ***'
86  if( io_l ) write(io_fid_log,*) '*** Data source : ', trim(h_source)
87  if( io_l ) write(io_fid_log,*) '*** Institute : ', trim(h_institute)
88 
89  allocate( axis_lon(imaxb,jmaxb) )
90  allocate( axis_lonx(imaxb,jmaxb) )
91  allocate( axis_lony(imaxb,jmaxb) )
92  allocate( axis_lonxy(imaxb,jmaxb) )
93  allocate( axis_lat(imaxb,jmaxb) )
94  allocate( axis_latx(imaxb,jmaxb) )
95  allocate( axis_laty(imaxb,jmaxb) )
96  allocate( axis_latxy(imaxb,jmaxb) )
97 
98  return
99  end subroutine fileio_setup
100 
101  !-----------------------------------------------------------------------------
103  subroutine fileio_set_coordinates( &
104  LON, &
105  LONX, &
106  LONY, &
107  LONXY, &
108  LAT, &
109  LATX, &
110  LATY, &
111  LATXY, &
112  CZ, &
113  FZ )
114  use scale_const, only: &
115  d2r => const_d2r
116  implicit none
117 
118  real(RP), intent(in) :: LON (ia,ja)
119  real(RP), intent(in) :: LONX(ia,ja)
120  real(RP), intent(in) :: LONY(ia,ja)
121  real(RP), intent(in) :: LONXY(ia,ja)
122  real(RP), intent(in) :: LAT (ia,ja)
123  real(RP), intent(in) :: LATX(ia,ja)
124  real(RP), intent(in) :: LATY(ia,ja)
125  real(RP), intent(in) :: LATXY(ia,ja)
126  real(RP), intent(in) :: CZ ( ka,ia,ja)
127  real(RP), intent(in) :: FZ (0:ka,ia,ja)
128  !---------------------------------------------------------------------------
129 
130  axis_lon(:,:) = lon(isb:ieb,jsb:jeb) / d2r
131  axis_lonx(:,:) = lonx(isb:ieb,jsb:jeb) / d2r
132  axis_lony(:,:) = lony(isb:ieb,jsb:jeb) / d2r
133  axis_lonxy(:,:) = lonxy(isb:ieb,jsb:jeb) / d2r
134  axis_lat(:,:) = lat(isb:ieb,jsb:jeb) / d2r
135  axis_latx(:,:) = latx(isb:ieb,jsb:jeb) / d2r
136  axis_laty(:,:) = laty(isb:ieb,jsb:jeb) / d2r
137  axis_latxy(:,:) = latxy(isb:ieb,jsb:jeb) / d2r
138 
139  return
140  end subroutine fileio_set_coordinates
141 
142  !-----------------------------------------------------------------------------
144  subroutine fileio_set_axes( &
145  fid, &
146  dtype, &
147  xy )
148  use gtool_file, only: &
149  fileputaxis, &
150  filesettattr, &
151  fileputassociatedcoordinates
152  use scale_grid, only: &
153  grid_cz, &
154  grid_cx, &
155  grid_cy, &
156  grid_fz, &
157  grid_fx, &
158  grid_fy, &
159  grid_cdz, &
160  grid_cdx, &
161  grid_cdy, &
162  grid_fdz, &
163  grid_fdx, &
164  grid_fdy, &
165  grid_cbfz, &
166  grid_cbfx, &
167  grid_cbfy, &
168  grid_fbfz, &
169  grid_fbfx, &
170  grid_fbfy, &
171  grid_cxg, &
172  grid_cyg, &
173  grid_fxg, &
174  grid_fyg, &
175  grid_cbfxg, &
176  grid_cbfyg, &
177  grid_fbfxg, &
178  grid_fbfyg
179  use scale_land_grid, only: &
180  grid_lcz, &
181  grid_lfz
182  use scale_urban_grid, only: &
183  grid_ucz, &
184  grid_ufz
185  implicit none
186 
187  integer, intent(in) :: fid
188  integer, intent(in) :: dtype
189  logical, intent(in), optional :: xy
190 
191  character(len=2) :: AXIS_name(2)
192  logical :: xy_
193  !---------------------------------------------------------------------------
194 
195  if ( present(xy) ) then
196  xy_ = xy
197  else
198  xy_ = .false.
199  end if
200 
201  if ( .NOT. xy_ ) then
202  call fileputaxis( fid, 'z', 'Z', 'm', 'z', dtype, grid_cz(ks:ke) )
203  end if
204  call fileputaxis( fid, 'x', 'X', 'm', 'x', dtype, grid_cx(isb:ieb) )
205  call fileputaxis( fid, 'y', 'Y', 'm', 'y', dtype, grid_cy(jsb:jeb) )
206  if ( .NOT. xy_ ) then
207  call fileputaxis( fid, 'zh', 'Z (half level)', 'm', 'zh', dtype, grid_fz(ks:ke) )
208  end if
209  call fileputaxis( fid, 'xh', 'X (half level)', 'm', 'xh', dtype, grid_fx(isb:ieb) )
210  call fileputaxis( fid, 'yh', 'Y (half level)', 'm', 'yh', dtype, grid_fy(jsb:jeb) )
211 
212  if ( .NOT. xy_ ) then
213  call fileputaxis( fid, 'lz', 'LZ', 'm', 'lz', dtype, grid_lcz(lks:lke) )
214  call fileputaxis( fid, 'lzh', 'LZ (half level)', 'm', 'lzh', dtype, grid_lfz(lks:lke) )
215  call fileputaxis( fid, 'uz', 'UZ', 'm', 'uz', dtype, grid_ucz(uks:uke) )
216  call fileputaxis( fid, 'uzh', 'UZ (half level)', 'm', 'uzh', dtype, grid_ufz(uks:uke) )
217  end if
218 
219 
220  if ( .NOT. xy_ ) then
221  call fileputaxis( fid, 'CZ', 'Atmos Grid Center Position Z', 'm', 'CZ', dtype, grid_cz )
222  end if
223  call fileputaxis( fid, 'CX', 'Atmos Grid Center Position X', 'm', 'CX', dtype, grid_cx )
224  call fileputaxis( fid, 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', dtype, grid_cy )
225  if ( .NOT. xy_ ) then
226  call fileputaxis( fid, 'FZ', 'Atmos Grid Face Position Z', 'm', 'FZ', dtype, grid_fz )
227  end if
228  call fileputaxis( fid, 'FX', 'Atmos Grid Face Position X', 'm', 'FX', dtype, grid_fx )
229  call fileputaxis( fid, 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', dtype, grid_fy )
230 
231  if ( .NOT. xy_ ) then
232  call fileputaxis( fid, 'CDZ', 'Grid Cell length Z', 'm', 'CZ', dtype, grid_cdz )
233  end if
234  call fileputaxis( fid, 'CDX', 'Grid Cell length X', 'm', 'CX', dtype, grid_cdx )
235  call fileputaxis( fid, 'CDY', 'Grid Cell length Y', 'm', 'CY', dtype, grid_cdy )
236  if ( .NOT. xy_ ) then
237  call fileputaxis( fid, 'FDZ', 'Grid distance Z', 'm', 'FDZ', dtype, grid_fdz )
238  end if
239  call fileputaxis( fid, 'FDX', 'Grid distance X', 'm', 'FDX', dtype, grid_fdx )
240  call fileputaxis( fid, 'FDY', 'Grid distance Y', 'm', 'FDY', dtype, grid_fdy )
241 
242  if ( .NOT. xy_ ) then
243  call fileputaxis( fid, 'LCZ', 'Land Grid Center Position Z', 'm', 'LCZ', dtype, grid_lcz )
244  call fileputaxis( fid, 'LFZ', 'Land Grid Face Position Z', 'm', 'LFZ', dtype, grid_lfz )
245  call fileputaxis( fid, 'LCDZ', 'Land Grid Cell length Z', 'm', 'LCZ', dtype, grid_lcz )
246 
247  call fileputaxis( fid, 'UCZ', 'Urban Grid Center Position Z', 'm', 'UCZ', dtype, grid_ucz )
248  call fileputaxis( fid, 'UFZ', 'Urban Grid Face Position Z', 'm', 'UFZ', dtype, grid_ufz )
249  call fileputaxis( fid, 'UCDZ', 'Urban Grid Cell length Z', 'm', 'UCZ', dtype, grid_ucz )
250  end if
251 
252  if ( .NOT. xy_ ) then
253  call fileputaxis( fid, 'CBFZ', 'Boundary factor Center Z', '1', 'CZ', dtype, grid_cbfz )
254  end if
255  call fileputaxis( fid, 'CBFX', 'Boundary factor Center X', '1', 'CX', dtype, grid_cbfx )
256  call fileputaxis( fid, 'CBFY', 'Boundary factor Center Y', '1', 'CY', dtype, grid_cbfy )
257  if ( .NOT. xy_ ) then
258  call fileputaxis( fid, 'FBFZ', 'Boundary factor Face Z', '1', 'CZ', dtype, grid_fbfz )
259  end if
260  call fileputaxis( fid, 'FBFX', 'Boundary factor Face X', '1', 'CX', dtype, grid_fbfx )
261  call fileputaxis( fid, 'FBFY', 'Boundary factor Face Y', '1', 'CY', dtype, grid_fbfy )
262 
263  call fileputaxis( fid, 'CXG', 'Grid Center Position X (global)', 'm', 'CXG', dtype, grid_cxg )
264  call fileputaxis( fid, 'CYG', 'Grid Center Position Y (global)', 'm', 'CYG', dtype, grid_cyg )
265  call fileputaxis( fid, 'FXG', 'Grid Face Position X (global)', 'm', 'FXG', dtype, grid_fxg )
266  call fileputaxis( fid, 'FYG', 'Grid Face Position Y (global)', 'm', 'FYG', dtype, grid_fyg )
267 
268  call fileputaxis( fid, 'CBFXG', 'Boundary factor Center X (global)', '1', 'CXG', dtype, grid_cbfxg )
269  call fileputaxis( fid, 'CBFYG', 'Boundary factor Center Y (global)', '1', 'CYG', dtype, grid_cbfyg )
270  call fileputaxis( fid, 'FBFXG', 'Boundary factor Face X (global)', '1', 'CXG', dtype, grid_fbfxg )
271  call fileputaxis( fid, 'FBFYG', 'Boundary factor Face Y (global)', '1', 'CYG', dtype, grid_fbfyg )
272 
273 
274  ! associate coordinates
275  axis_name = (/'x ','y '/)
276  call fileputassociatedcoordinates( fid, 'lon' , 'longitude' , &
277  'degrees_east' , axis_name, dtype, axis_lon(:,:) )
278  axis_name = (/'xh','y '/)
279  call fileputassociatedcoordinates( fid, 'lon_uy', 'longitude (half level uy)', &
280  'degrees_east' , axis_name, dtype, axis_lonx(:,:) )
281  axis_name = (/'x ','yh'/)
282  call fileputassociatedcoordinates( fid, 'lon_xv', 'longitude (half level xv)', &
283  'degrees_east' , axis_name, dtype, axis_lony(:,:) )
284  axis_name = (/'xh','yh'/)
285  call fileputassociatedcoordinates( fid, 'lon_uv', 'longitude (half level uv)', &
286  'degrees_east' , axis_name, dtype, axis_lonxy(:,:) )
287  axis_name = (/'x ','y '/)
288  call fileputassociatedcoordinates( fid, 'lat' , 'latitude' , &
289  'degrees_north', axis_name, dtype, axis_lat(:,:) )
290  axis_name = (/'xh','y '/)
291  call fileputassociatedcoordinates( fid, 'lat_uy', 'latitude (half level uy)' , &
292  'degrees_north', axis_name, dtype, axis_latx(:,:) )
293  axis_name = (/'x ','yh'/)
294  call fileputassociatedcoordinates( fid, 'lat_xv', 'latitude (half level xv)' , &
295  'degrees_north', axis_name, dtype, axis_laty(:,:) )
296  axis_name = (/'xh','yh'/)
297  call fileputassociatedcoordinates( fid, 'lat_uv', 'latitude (half level uv)' , &
298  'degrees_north', axis_name, dtype, axis_latxy(:,:) )
299 
300  ! attributes
301  if ( .NOT. xy_ ) then
302  call filesettattr( fid, 'lz', 'positive', 'down' )
303  call filesettattr( fid, 'lzh', 'positive', 'down' )
304  call filesettattr( fid, 'uz', 'positive', 'down' )
305  call filesettattr( fid, 'uzh', 'positive', 'down' )
306  call filesettattr( fid, 'LCZ', 'positive', 'down' )
307  call filesettattr( fid, 'LFZ', 'positive', 'down' )
308  call filesettattr( fid, 'UCZ', 'positive', 'down' )
309  call filesettattr( fid, 'UFZ', 'positive', 'down' )
310  end if
311 
312  return
313  end subroutine fileio_set_axes
314 
315  !-----------------------------------------------------------------------------
317  subroutine fileio_read_1d( &
318  var, &
319  basename, &
320  varname, &
321  axistype, &
322  step )
323  use gtool_file, only: &
324  fileread
325  use scale_process, only: &
326  prc_myrank, &
328  implicit none
329 
330  real(RP), intent(out) :: var(:)
331  character(len=*), intent(in) :: basename
332  character(len=*), intent(in) :: varname
333  character(len=*), intent(in) :: axistype
334  integer, intent(in) :: step
335 
336  integer :: dim1_max, dim1_S, dim1_E
337  real(RP), allocatable :: var1D(:)
338  !---------------------------------------------------------------------------
339 
340  call prof_rapstart('FILE_I_NetCDF', 2)
341 
342  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 1D var: ', trim(varname)
343 
344  if ( axistype == 'Z' ) then
345  dim1_max = kmax
346  dim1_s = ks
347  dim1_e = ke
348  elseif( axistype == 'X' ) then
349  dim1_max = imaxb
350  dim1_s = isb
351  dim1_e = ieb
352  elseif( axistype == 'Y' ) then
353  dim1_max = jmaxb
354  dim1_s = jsb
355  dim1_e = jeb
356  else
357  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
358  call prc_mpistop
359  endif
360 
361  allocate( var1d(dim1_max) )
362 
363  call fileread( var1d(:), basename, varname, step, prc_myrank )
364  var(dim1_s:dim1_e) = var1d(1:dim1_max)
365 
366  deallocate( var1d )
367 
368  call prof_rapend ('FILE_I_NetCDF', 2)
369 
370  return
371  end subroutine fileio_read_1d
372 
373  !-----------------------------------------------------------------------------
375  subroutine fileio_read_2d( &
376  var, &
377  basename, &
378  varname, &
379  axistype, &
380  step )
381  use gtool_file, only: &
382  fileread
383  use scale_process, only: &
384  prc_myrank, &
386  implicit none
387 
388  real(RP), intent(out) :: var(:,:)
389  character(len=*), intent(in) :: basename
390  character(len=*), intent(in) :: varname
391  character(len=*), intent(in) :: axistype
392  integer, intent(in) :: step
393 
394  integer :: dim1_max, dim1_S, dim1_E
395  integer :: dim2_max, dim2_S, dim2_E
396  real(RP), allocatable :: var2D(:,:)
397  !---------------------------------------------------------------------------
398 
399  call prof_rapstart('FILE_I_NetCDF', 2)
400 
401  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 2D var: ', trim(varname)
402 
403  if ( axistype == 'XY' ) then
404  dim1_max = imaxb
405  dim2_max = jmaxb
406  dim1_s = isb
407  dim1_e = ieb
408  dim2_s = jsb
409  dim2_e = jeb
410  elseif( axistype == 'ZX' ) then
411  dim1_max = kmax
412  dim2_max = imaxb
413  dim1_s = ks
414  dim1_e = ke
415  dim2_s = isb
416  dim2_e = ieb
417  else
418  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
419  call prc_mpistop
420  endif
421 
422  allocate( var2d(dim1_max,dim2_max) )
423 
424  call fileread( var2d(:,:), basename, varname, step, prc_myrank )
425  var(dim1_s:dim1_e,dim2_s:dim2_e) = var2d(1:dim1_max,1:dim2_max)
426 
427  deallocate( var2d )
428 
429  call prof_rapend ('FILE_I_NetCDF', 2)
430 
431  return
432  end subroutine fileio_read_2d
433 
434  !-----------------------------------------------------------------------------
436  subroutine fileio_read_3d( &
437  var, &
438  basename, &
439  varname, &
440  axistype, &
441  step )
442  use gtool_file, only: &
443  fileread
444  use scale_process, only: &
445  prc_myrank, &
447  implicit none
448 
449  real(RP), intent(out) :: var(:,:,:)
450  character(len=*), intent(in) :: basename
451  character(len=*), intent(in) :: varname
452  character(len=*), intent(in) :: axistype
453  integer, intent(in) :: step
454 
455  integer :: dim1_max, dim1_S, dim1_E
456  integer :: dim2_max, dim2_S, dim2_E
457  integer :: dim3_max, dim3_S, dim3_E
458  real(RP), allocatable :: var3D(:,:,:)
459  !---------------------------------------------------------------------------
460 
461  call prof_rapstart('FILE_I_NetCDF', 2)
462 
463  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 3D var: ', trim(varname)
464 
465  if ( axistype == 'ZXY' ) then
466  dim1_max = kmax
467  dim2_max = imaxb
468  dim3_max = jmaxb
469  dim1_s = ks
470  dim1_e = ke
471  dim2_s = isb
472  dim2_e = ieb
473  dim3_s = jsb
474  dim3_e = jeb
475  else if ( axistype == 'XYT' ) then
476  dim1_max = imaxb
477  dim2_max = jmaxb
478  dim3_max = step
479  dim1_s = isb
480  dim1_e = ieb
481  dim2_s = jsb
482  dim2_e = jeb
483  dim3_s = 1
484  dim3_e = step
485  else if ( axistype == 'Land' ) then
486  dim1_max = lkmax
487  dim2_max = imaxb
488  dim3_max = jmaxb
489  dim1_s = lks
490  dim1_e = lke
491  dim2_s = isb
492  dim2_e = ieb
493  dim3_s = jsb
494  dim3_e = jeb
495  else if ( axistype == 'Urban' ) then
496  dim1_max = ukmax
497  dim2_max = imaxb
498  dim3_max = jmaxb
499  dim1_s = uks
500  dim1_e = uke
501  dim2_s = isb
502  dim2_e = ieb
503  dim3_s = jsb
504  dim3_e = jeb
505  else
506  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
507  call prc_mpistop
508  endif
509 
510  allocate( var3d(dim1_max,dim2_max,dim3_max) )
511 
512  call fileread( var3d(:,:,:), basename, varname, step, prc_myrank )
513  var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e) = var3d(1:dim1_max,1:dim2_max,1:dim3_max)
514 
515  deallocate( var3d )
516 
517  call prof_rapend ('FILE_I_NetCDF', 2)
518 
519  return
520  end subroutine fileio_read_3d
521 
522  !-----------------------------------------------------------------------------
524  subroutine fileio_read_4d( &
525  var, &
526  basename, &
527  varname, &
528  axistype, &
529  step )
530  use gtool_file, only: &
531  fileread
532  use scale_process, only: &
533  prc_myrank, &
535  implicit none
536 
537  real(RP), intent(out) :: var(:,:,:,:)
538  character(len=*), intent(in) :: basename
539  character(len=*), intent(in) :: varname
540  character(len=*), intent(in) :: axistype
541  integer, intent(in) :: step
542 
543  integer :: dim1_max, dim1_S, dim1_E
544  integer :: dim2_max, dim2_S, dim2_E
545  integer :: dim3_max, dim3_S, dim3_E
546  integer :: dim4_max, dim4_S, dim4_E
547  real(RP), allocatable :: var4D(:,:,:,:)
548  !---------------------------------------------------------------------------
549 
550  call prof_rapstart('FILE_I_NetCDF', 2)
551 
552  if( io_l ) write(io_fid_log,'(1x,A,A15)') '*** Read 4D var: ', trim(varname)
553 
554  if ( axistype == 'ZXYT' ) then
555  dim1_max = kmax
556  dim2_max = imaxb
557  dim3_max = jmaxb
558  dim4_max = step
559  dim1_s = ks
560  dim1_e = ke
561  dim2_s = isb
562  dim2_e = ieb
563  dim3_s = jsb
564  dim3_e = jeb
565  dim4_s = 1
566  dim4_e = step
567  else
568  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
569  call prc_mpistop
570  endif
571 
572  allocate( var4d(dim1_max,dim2_max,dim3_max,dim4_max) )
573 
574  call fileread( var4d(:,:,:,:), basename, varname, step, prc_myrank )
575  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)
576 
577  deallocate( var4d )
578 
579  call prof_rapend ('FILE_I_NetCDF', 2)
580 
581  return
582  end subroutine fileio_read_4d
583 
584  !-----------------------------------------------------------------------------
586  subroutine fileio_write_1d( &
587  var, &
588  basename, &
589  title, &
590  varname, &
591  desc, &
592  unit, &
593  axistype, &
594  datatype, &
595  date, &
596  subsec, &
597  append )
598  use gtool_file_h, only: &
599  file_real8, &
600  file_real4
601  use gtool_file, only: &
602  filecreate, &
603  fileaddvariable, &
604  filesetglobalattribute, &
605  filewrite
606  use scale_process, only: &
607  prc_masterrank, &
608  prc_myrank, &
610  use scale_rm_process, only: &
611  prc_2drank
612  use scale_time, only: &
613  nowdate => time_nowdate, &
614  nowms => time_nowms, &
615  nowsec => time_nowdaysec
616  implicit none
617 
618  real(RP), intent(in) :: var(:)
619  character(len=*), intent(in) :: basename
620  character(len=*), intent(in) :: title
621  character(len=*), intent(in) :: varname
622  character(len=*), intent(in) :: desc
623  character(len=*), intent(in) :: unit
624  character(len=*), intent(in) :: axistype
625  character(len=*), intent(in) :: datatype
626 
627  integer, optional, intent(in) :: date(6)
628  real(DP),optional, intent(in) :: subsec
629  logical, optional, intent(in) :: append
630 
631  integer :: dtype
632  character(len=2) :: dims(1)
633  integer :: dim1_max, dim1_S, dim1_E
634  real(RP), allocatable :: var1D(:)
635 
636  integer :: rankidx(2)
637  logical :: fileexisted
638  integer :: fid, vid
639  character(len=34) :: tunits
640  !---------------------------------------------------------------------------
641 
642  call prof_rapstart('FILE_O_NetCDF', 2)
643 
644  rankidx(1) = prc_2drank(prc_myrank,1)
645  rankidx(2) = prc_2drank(prc_myrank,2)
646 
647  if ( datatype == 'REAL8' ) then
648  dtype = file_real8
649  elseif( datatype == 'REAL4' ) then
650  dtype = file_real4
651  else
652  if ( rp == 8 ) then
653  dtype = file_real8
654  elseif( rp == 4 ) then
655  dtype = file_real4
656  else
657  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
658  call prc_mpistop
659  endif
660  endif
661 
662  call filecreate( fid, & ! [OUT]
663  fileexisted, & ! [OUT]
664  basename, & ! [IN]
665  title, & ! [IN]
666  h_source, & ! [IN]
667  h_institute, & ! [IN]
668  prc_masterrank, & ! [IN]
669  prc_myrank, & ! [IN]
670  rankidx, & ! [IN]
671  append = append ) ! [IN]
672 
673  if ( .NOT. fileexisted ) then ! only once
674  call fileio_set_axes( fid, dtype ) ! [IN]
675  if ( present( subsec ) ) then
676  call filesetglobalattribute( fid, "time", (/subsec/) )
677  else
678  call filesetglobalattribute( fid, "time", (/nowms/) )
679  end if
680  if ( present( date ) ) then
681  call getcftunits(tunits, date)
682  else
683  call getcftunits(tunits, nowdate)
684  end if
685  call filesetglobalattribute( fid, "time_units", tunits )
686  endif
687 
688  if ( axistype == 'Z' ) then
689  dims = (/'z'/)
690  dim1_max = kmax
691  dim1_s = ks
692  dim1_e = ke
693  elseif( axistype == 'X' ) then
694  dims = (/'x'/)
695  dim1_max = imaxb
696  dim1_s = isb
697  dim1_e = ieb
698  elseif( axistype == 'Y' ) then
699  dims = (/'y'/)
700  dim1_max = jmaxb
701  dim1_s = jsb
702  dim1_e = jeb
703  else
704  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
705  call prc_mpistop
706  endif
707 
708  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype ) ! [IN]
709 
710  allocate( var1d(dim1_max) )
711 
712  var1d(1:dim1_max) = var(dim1_s:dim1_e)
713  call filewrite( fid, vid, var1d(:), nowsec, nowsec ) ! [IN]
714 
715  deallocate( var1d )
716 
717  call prof_rapend ('FILE_O_NetCDF', 2)
718 
719  return
720  end subroutine fileio_write_1d
721 
722  !-----------------------------------------------------------------------------
724  subroutine fileio_write_2d( &
725  var, &
726  basename, &
727  title, &
728  varname, &
729  desc, &
730  unit, &
731  axistype, &
732  datatype, &
733  date, &
734  subsec, &
735  append, &
736  nohalo, &
737  nozcoord )
738  use gtool_file, only: &
739  rmiss
740  use gtool_file_h, only: &
741  file_real8, &
742  file_real4
743  use gtool_file, only: &
744  filecreate, &
745  fileaddvariable, &
746  filesetglobalattribute, &
747  filewrite
748  use scale_process, only: &
749  prc_masterrank, &
750  prc_myrank, &
752  use scale_rm_process, only: &
753  prc_2drank
754  use scale_time, only: &
755  nowdate => time_nowdate, &
756  nowms => time_nowms, &
757  nowsec => time_nowdaysec
758  implicit none
759 
760  real(RP), intent(in) :: var(:,:)
761  character(len=*), intent(in) :: basename
762  character(len=*), intent(in) :: title
763  character(len=*), intent(in) :: varname
764  character(len=*), intent(in) :: desc
765  character(len=*), intent(in) :: unit
766  character(len=*), intent(in) :: axistype
767  character(len=*), intent(in) :: datatype
768  integer, optional, intent(in) :: date(6)
769  real(DP),optional, intent(in) :: subsec
770  logical, optional, intent(in) :: append
771  logical, optional, intent(in) :: nohalo
772  logical, optional, intent(in) :: nozcoord
773 
774  real(RP) :: varhalo( size(var(:,1)), size(var(1,:)) )
775 
776  integer :: dtype
777  character(len=2) :: dims(2)
778  integer :: dim1_max, dim1_S, dim1_E
779  integer :: dim2_max, dim2_S, dim2_E
780  real(RP), allocatable :: var2D(:,:)
781 
782  integer :: rankidx(2)
783  logical :: fileexisted
784  integer :: fid, vid
785  integer :: i, j
786  logical :: nohalo_
787  character(len=34) :: tunits
788  !---------------------------------------------------------------------------
789 
790  call prof_rapstart('FILE_O_NetCDF', 2)
791 
792  nohalo_ = .false.
793  if ( present(nohalo) ) nohalo_ = nohalo
794 
795  rankidx(1) = prc_2drank(prc_myrank,1)
796  rankidx(2) = prc_2drank(prc_myrank,2)
797 
798  if ( datatype == 'REAL8' ) then
799  dtype = file_real8
800  elseif( datatype == 'REAL4' ) then
801  dtype = file_real4
802  else
803  if ( rp == 8 ) then
804  dtype = file_real8
805  elseif( rp == 4 ) then
806  dtype = file_real4
807  else
808  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
809  call prc_mpistop
810  endif
811  endif
812 
813  call filecreate( fid, & ! [OUT]
814  fileexisted, & ! [OUT]
815  basename, & ! [IN]
816  title, & ! [IN]
817  h_source, & ! [IN]
818  h_institute, & ! [IN]
819  prc_masterrank, & ! [IN]
820  prc_myrank, & ! [IN]
821  rankidx, & ! [IN]
822  append = append ) ! [IN]
823 
824  if ( axistype == 'XY' ) then
825  dims = (/'x','y'/)
826  dim1_max = imaxb
827  dim2_max = jmaxb
828  dim1_s = isb
829  dim1_e = ieb
830  dim2_s = jsb
831  dim2_e = jeb
832  elseif ( axistype == 'UY' ) then
833  dims = (/'xh','y '/)
834  dim1_max = imaxb
835  dim2_max = jmaxb
836  dim1_s = isb
837  dim1_e = ieb
838  dim2_s = jsb
839  dim2_e = jeb
840  elseif ( axistype == 'XV' ) then
841  dims = (/'x ','yh'/)
842  dim1_max = imaxb
843  dim2_max = jmaxb
844  dim1_s = isb
845  dim1_e = ieb
846  dim2_s = jsb
847  dim2_e = jeb
848  elseif ( axistype == 'UV' ) then
849  dims = (/'xh','yh'/)
850  dim1_max = imaxb
851  dim2_max = jmaxb
852  dim1_s = isb
853  dim1_e = ieb
854  dim2_s = jsb
855  dim2_e = jeb
856  elseif( axistype == 'ZX' ) then
857  dims = (/'z','x'/)
858  dim1_max = kmax
859  dim2_max = imaxb
860  dim1_s = ks
861  dim1_e = ke
862  dim2_s = isb
863  dim2_e = ieb
864  else
865  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
866  call prc_mpistop
867  endif
868 
869  if ( .NOT. fileexisted ) then ! only once
870  call fileio_set_axes( fid, dtype, nozcoord ) ! [IN]
871  if ( present( subsec ) ) then
872  call filesetglobalattribute( fid, "time", (/subsec/) )
873  else
874  call filesetglobalattribute( fid, "time", (/nowms/) )
875  end if
876  if ( present( date ) ) then
877  call getcftunits(tunits, date)
878  else
879  call getcftunits(tunits, nowdate)
880  end if
881  call filesetglobalattribute( fid, "time_units", tunits )
882  endif
883 
884  varhalo(:,:) = var(:,:)
885 
886  if ( nohalo_ ) then
887  ! W halo
888  do j = 1, ja
889  do i = 1, is-1
890  varhalo(i,j) = rmiss
891  end do
892  end do
893  ! E halo
894  do j = 1, ja
895  do i = ie+1, ia
896  varhalo(i,j) = rmiss
897  end do
898  end do
899  ! S halo
900  do j = 1, js-1
901  do i = 1, ia
902  varhalo(i,j) = rmiss
903  end do
904  end do
905  ! N halo
906  do j = je+1, ja
907  do i = 1, ia
908  varhalo(i,j) = rmiss
909  end do
910  end do
911  end if
912 
913  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype ) ! [IN]
914 
915  allocate( var2d(dim1_max,dim2_max) )
916 
917  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
918  call filewrite( fid, vid, var2d(:,:), nowsec, nowsec ) ! [IN]
919 
920  deallocate( var2d )
921 
922  call prof_rapend ('FILE_O_NetCDF', 2)
923 
924  return
925  end subroutine fileio_write_2d
926 
927  !-----------------------------------------------------------------------------
929  subroutine fileio_write_3d( &
930  var, &
931  basename, &
932  title, &
933  varname, &
934  desc, &
935  unit, &
936  axistype, &
937  datatype, &
938  date, &
939  subsec, &
940  append, &
941  nohalo )
942  use gtool_file, only: &
943  rmiss
944  use gtool_file_h, only: &
945  file_real8, &
946  file_real4
947  use gtool_file, only: &
948  filecreate, &
949  fileaddvariable, &
950  filesetglobalattribute, &
951  filewrite
952  use scale_process, only: &
953  prc_masterrank, &
954  prc_myrank, &
956  use scale_rm_process, only: &
957  prc_2drank
958  use scale_time, only: &
959  nowdate => time_nowdate, &
960  nowms => time_nowms, &
961  nowsec => time_nowdaysec
962  implicit none
963 
964  real(RP), intent(in) :: var(:,:,:)
965  character(len=*), intent(in) :: basename
966  character(len=*), intent(in) :: title
967  character(len=*), intent(in) :: varname
968  character(len=*), intent(in) :: desc
969  character(len=*), intent(in) :: unit
970  character(len=*), intent(in) :: axistype
971  character(len=*), intent(in) :: datatype
972  integer, optional, intent(in) :: date(6)
973  real(DP),optional, intent(in) :: subsec
974  logical, optional, intent(in) :: append
975  logical, optional, intent(in) :: nohalo
976 
977  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)), size(var(1,1,:)) )
978 
979  integer :: dtype
980  character(len=2) :: dims(3)
981  integer :: dim1_max, dim1_S, dim1_E
982  integer :: dim2_max, dim2_S, dim2_E
983  integer :: dim3_max, dim3_S, dim3_E
984 
985  real(RP), allocatable :: var3D(:,:,:)
986 
987  integer :: rankidx(2)
988  logical :: append_sw
989  logical :: fileexisted
990  integer :: fid, vid
991  integer :: i, j, k
992  logical :: nohalo_
993  character(len=34) :: tunits
994  !---------------------------------------------------------------------------
995 
996  call prof_rapstart('FILE_O_NetCDF', 2)
997 
998  nohalo_ = .false.
999  if ( present(nohalo) ) nohalo_ = nohalo
1000 
1001  rankidx(1) = prc_2drank(prc_myrank,1)
1002  rankidx(2) = prc_2drank(prc_myrank,2)
1003 
1004  if ( datatype == 'REAL8' ) then
1005  dtype = file_real8
1006  elseif( datatype == 'REAL4' ) then
1007  dtype = file_real4
1008  else
1009  if ( rp == 8 ) then
1010  dtype = file_real8
1011  elseif( rp == 4 ) then
1012  dtype = file_real4
1013  else
1014  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
1015  call prc_mpistop
1016  endif
1017  endif
1018 
1019  append_sw = .false.
1020  if ( present(append) ) then
1021  append_sw = append
1022  endif
1023 
1024  call filecreate( fid, & ! [OUT]
1025  fileexisted, & ! [OUT]
1026  basename, & ! [IN]
1027  title, & ! [IN]
1028  h_source, & ! [IN]
1029  h_institute, & ! [IN]
1030  prc_masterrank, & ! [IN]
1031  prc_myrank, & ! [IN]
1032  rankidx, & ! [IN]
1033  append = append_sw ) ! [IN]
1034 
1035  if ( .NOT. fileexisted ) then ! only once
1036  call fileio_set_axes( fid, dtype ) ! [IN]
1037  if ( present( subsec ) ) then
1038  call filesetglobalattribute( fid, "time", (/subsec/) )
1039  else
1040  call filesetglobalattribute( fid, "time", (/nowms/) )
1041  end if
1042  if ( present( date ) ) then
1043  call getcftunits(tunits, date)
1044  else
1045  call getcftunits(tunits, nowdate)
1046  end if
1047  call filesetglobalattribute( fid, "time_units", tunits )
1048  endif
1049 
1050  if ( axistype == 'ZXY' ) then
1051  dims = (/'z','x','y'/)
1052  dim1_max = kmax
1053  dim2_max = imaxb
1054  dim3_max = jmaxb
1055  dim1_s = ks
1056  dim1_e = ke
1057  dim2_s = isb
1058  dim2_e = ieb
1059  dim3_s = jsb
1060  dim3_e = jeb
1061  elseif( axistype == 'ZHXY' ) then
1062  dims = (/'zh','x ','y '/)
1063  dim1_max = kmax
1064  dim2_max = imaxb
1065  dim3_max = jmaxb
1066  dim1_s = ks
1067  dim1_e = ke
1068  dim2_s = isb
1069  dim2_e = ieb
1070  dim3_s = jsb
1071  dim3_e = jeb
1072  elseif( axistype == 'ZXHY' ) then
1073  dims = (/'z ','xh','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 == 'ZXYH' ) then
1084  dims = (/'z ','x ','yh'/)
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 == 'Land' ) then
1095  dims = (/'lz','x ','y '/)
1096  dim1_max = lkmax
1097  dim2_max = imaxb
1098  dim3_max = jmaxb
1099  dim1_s = lks
1100  dim1_e = lke
1101  dim2_s = isb
1102  dim2_e = ieb
1103  dim3_s = jsb
1104  dim3_e = jeb
1105  elseif( axistype == 'Urban' ) then
1106  dims = (/'uz','x ','y '/)
1107  dim1_max = ukmax
1108  dim2_max = imaxb
1109  dim3_max = jmaxb
1110  dim1_s = uks
1111  dim1_e = uke
1112  dim2_s = isb
1113  dim2_e = ieb
1114  dim3_s = jsb
1115  dim3_e = jeb
1116  else
1117  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
1118  call prc_mpistop
1119  endif
1120 
1121  varhalo(:,:,:) = var(:,:,:)
1122 
1123  if ( nohalo_ ) then
1124  ! W halo
1125  do k = 1, dim1_max
1126  do j = 1, ja
1127  do i = 1, is-1
1128  varhalo(k,i,j) = rmiss
1129  end do
1130  end do
1131  end do
1132  ! E halo
1133  do k = 1, dim1_max
1134  do j = 1, ja
1135  do i = ie+1, ia
1136  varhalo(k,i,j) = rmiss
1137  end do
1138  end do
1139  end do
1140  ! S halo
1141  do k = 1, dim1_max
1142  do j = 1, js-1
1143  do i = 1, ia
1144  varhalo(k,i,j) = rmiss
1145  end do
1146  end do
1147  end do
1148  ! N halo
1149  do k = 1, dim1_max
1150  do j = je+1, ja
1151  do i = 1, ia
1152  varhalo(k,i,j) = rmiss
1153  end do
1154  end do
1155  end do
1156  end if
1157 
1158  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype ) ! [IN]
1159 
1160  allocate( var3d(dim1_max,dim2_max,dim3_max) )
1161 
1162  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
1163  call filewrite( fid, vid, var3d(:,:,:), nowsec, nowsec ) ! [IN]
1164 
1165  deallocate( var3d )
1166 
1167  call prof_rapend ('FILE_O_NetCDF', 2)
1168 
1169  return
1170  end subroutine fileio_write_3d
1171 
1172  !-----------------------------------------------------------------------------
1174  subroutine fileio_write_3d_t( &
1175  var, &
1176  basename, &
1177  title, &
1178  varname, &
1179  desc, &
1180  unit, &
1181  axistype, &
1182  datatype, &
1183  timeintv, &
1184  tsince, &
1185  append, &
1186  timetarg, &
1187  nohalo )
1188  use gtool_file, only: &
1189  rmiss
1190  use gtool_file_h, only: &
1191  file_real8, &
1192  file_real4
1193  use gtool_file, only: &
1194  filecreate, &
1195  fileputaxis, &
1196  fileaddvariable, &
1197  filewrite
1198  use scale_process, only: &
1199  prc_masterrank, &
1200  prc_myrank, &
1201  prc_mpistop
1202  use scale_rm_process, only: &
1203  prc_2drank
1204  implicit none
1205 
1206  real(RP), intent(in) :: var(:,:,:)
1207  character(len=*), intent(in) :: basename
1208  character(len=*), intent(in) :: title
1209  character(len=*), intent(in) :: varname
1210  character(len=*), intent(in) :: desc
1211  character(len=*), intent(in) :: unit
1212  character(len=*), intent(in) :: axistype
1213  character(len=*), intent(in) :: datatype
1214  real(RP), intent(in) :: timeintv
1215  integer , intent(in) :: tsince(6)
1216  logical, optional, intent(in) :: append
1217  integer, optional, intent(in) :: timetarg
1218  logical, optional, intent(in) :: nohalo
1219 
1220  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)) )
1221 
1222  integer :: dtype
1223  character(len=2) :: dims(2)
1224  integer :: dim1_max, dim1_S, dim1_E
1225  integer :: dim2_max, dim2_S, dim2_E
1226 
1227  real(RP), allocatable :: var2D(:,:)
1228  real(DP) :: time_interval, nowtime
1229 
1230  character(len=34) :: tunits
1231 
1232  integer :: rankidx(2)
1233  logical :: append_sw
1234  logical :: fileexisted
1235  integer :: fid, vid
1236  integer :: step
1237  integer :: i, j, n
1238  logical :: nohalo_
1239  !---------------------------------------------------------------------------
1240 
1241  call prof_rapstart('FILE_O_NetCDF', 2)
1242 
1243  nohalo_ = .false.
1244  if ( present(nohalo) ) nohalo_ = nohalo
1245 
1246  time_interval = timeintv
1247  step = size(var(isb,jsb,:))
1248 
1249  rankidx(1) = prc_2drank(prc_myrank,1)
1250  rankidx(2) = prc_2drank(prc_myrank,2)
1251 
1252  if ( datatype == 'REAL8' ) then
1253  dtype = file_real8
1254  elseif( datatype == 'REAL4' ) then
1255  dtype = file_real4
1256  else
1257  if ( rp == 8 ) then
1258  dtype = file_real8
1259  elseif( rp == 4 ) then
1260  dtype = file_real4
1261  else
1262  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
1263  call prc_mpistop
1264  endif
1265  endif
1266 
1267  append_sw = .false.
1268  if ( present(append) ) then
1269  append_sw = append
1270  endif
1271 
1272  write(tunits,'(a,i4.4,"-",i2.2,"-",i2.2," ",i2.2,":",i2.2,":",i2.2)') 'seconds since ', tsince
1273 
1274  call filecreate( fid, & ! [OUT]
1275  fileexisted, & ! [OUT]
1276  basename, & ! [IN]
1277  title, & ! [IN]
1278  h_source, & ! [IN]
1279  h_institute, & ! [IN]
1280  prc_masterrank, & ! [IN]
1281  prc_myrank, & ! [IN]
1282  rankidx, & ! [IN]
1283  time_units = tunits, & ! [IN]
1284  append = append_sw ) ! [IN]
1285 
1286  if ( .NOT. fileexisted ) then ! only once
1287  call fileio_set_axes( fid, dtype ) ! [IN]
1288  endif
1289 
1290  if ( axistype == 'XYT' ) then
1291  dims = (/'x','y'/)
1292  dim1_max = imaxb
1293  dim2_max = jmaxb
1294  dim1_s = isb
1295  dim1_e = ieb
1296  dim2_s = jsb
1297  dim2_e = jeb
1298  else
1299  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
1300  call prc_mpistop
1301  endif
1302 
1303  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype, time_interval ) ! [IN]
1304  allocate( var2d(dim1_max,dim2_max) )
1305 
1306  if ( present(timetarg) ) then
1307  varhalo(:,:) = var(:,:,timetarg)
1308 
1309  if ( nohalo_ ) then
1310  ! W halo
1311  do j = 1, ja
1312  do i = 1, is-1
1313  varhalo(i,j) = rmiss
1314  end do
1315  end do
1316  ! E halo
1317  do j = 1, ja
1318  do i = ie+1, ia
1319  varhalo(i,j) = rmiss
1320  end do
1321  end do
1322  ! S halo
1323  do j = 1, js-1
1324  do i = 1, ia
1325  varhalo(i,j) = rmiss
1326  end do
1327  end do
1328  ! N halo
1329  do j = je+1, ja
1330  do i = 1, ia
1331  varhalo(i,j) = rmiss
1332  end do
1333  end do
1334  end if
1335 
1336  nowtime = (timetarg-1) * time_interval
1337  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
1338  call filewrite( fid, vid, var2d(:,:), nowtime, nowtime ) ! [IN]
1339  else
1340  nowtime = 0.0_dp
1341  do n = 1, step
1342  varhalo(:,:) = var(:,:,n)
1343 
1344  if ( nohalo_ ) then
1345  ! W halo
1346  do j = 1, ja
1347  do i = 1, is-1
1348  varhalo(i,j) = rmiss
1349  end do
1350  end do
1351  ! E halo
1352  do j = 1, ja
1353  do i = ie+1, ia
1354  varhalo(i,j) = rmiss
1355  end do
1356  end do
1357  ! S halo
1358  do j = 1, js-1
1359  do i = 1, ia
1360  varhalo(i,j) = rmiss
1361  end do
1362  end do
1363  ! N halo
1364  do j = je+1, ja
1365  do i = 1, ia
1366  varhalo(i,j) = rmiss
1367  end do
1368  end do
1369  end if
1370 
1371  var2d(1:dim1_max,1:dim2_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e)
1372  call filewrite( fid, vid, var2d(:,:), nowtime, nowtime ) ! [IN]
1373  nowtime = nowtime + time_interval
1374  enddo
1375  endif
1376 
1377  deallocate( var2d )
1378 
1379  call prof_rapend ('FILE_O_NetCDF', 2)
1380 
1381  return
1382  end subroutine fileio_write_3d_t
1383 
1384  !-----------------------------------------------------------------------------
1386  subroutine fileio_write_4d( &
1387  var, &
1388  basename, &
1389  title, &
1390  varname, &
1391  desc, &
1392  unit, &
1393  axistype, &
1394  datatype, &
1395  timeintv, &
1396  tsince, &
1397  append, &
1398  timetarg, &
1399  nohalo )
1400  use gtool_file, only: &
1401  rmiss
1402  use gtool_file_h, only: &
1403  file_real8, &
1404  file_real4
1405  use gtool_file, only: &
1406  filecreate, &
1407  fileputaxis, &
1408  fileaddvariable, &
1409  filewrite
1410  use scale_process, only: &
1411  prc_masterrank, &
1412  prc_myrank, &
1413  prc_mpistop
1414  use scale_rm_process, only: &
1415  prc_2drank
1416  implicit none
1417 
1418  real(RP), intent(in) :: var(:,:,:,:)
1419  character(len=*), intent(in) :: basename
1420  character(len=*), intent(in) :: title
1421  character(len=*), intent(in) :: varname
1422  character(len=*), intent(in) :: desc
1423  character(len=*), intent(in) :: unit
1424  character(len=*), intent(in) :: axistype
1425  character(len=*), intent(in) :: datatype
1426  real(RP), intent(in) :: timeintv
1427  integer , intent(in) :: tsince(6)
1428  logical, optional, intent(in) :: append
1429  integer, optional, intent(in) :: timetarg
1430  logical, optional, intent(in) :: nohalo
1431 
1432  real(RP) :: varhalo( size(var(:,1,1,1)), size(var(1,:,1,1)), size(var(1,1,:,1)) )
1433 
1434  integer :: dtype
1435  character(len=2) :: dims(3)
1436  integer :: dim1_max, dim1_S, dim1_E
1437  integer :: dim2_max, dim2_S, dim2_E
1438  integer :: dim3_max, dim3_S, dim3_E
1439 
1440  real(RP), allocatable :: var3D(:,:,:)
1441  real(DP) :: time_interval, nowtime
1442 
1443  character(len=34) :: tunits
1444 
1445  integer :: rankidx(2)
1446  logical :: append_sw
1447  logical :: fileexisted
1448  integer :: fid, vid
1449  integer :: step
1450  integer :: i, j, k, n
1451  logical :: nohalo_
1452  !---------------------------------------------------------------------------
1453 
1454  call prof_rapstart('FILE_O_NetCDF', 2)
1455 
1456  nohalo_ = .false.
1457  if ( present(nohalo) ) nohalo_ = nohalo
1458 
1459  time_interval = timeintv
1460  step = size(var(ks,isb,jsb,:))
1461 
1462  rankidx(1) = prc_2drank(prc_myrank,1)
1463  rankidx(2) = prc_2drank(prc_myrank,2)
1464 
1465  if ( datatype == 'REAL8' ) then
1466  dtype = file_real8
1467  elseif( datatype == 'REAL4' ) then
1468  dtype = file_real4
1469  else
1470  if ( rp == 8 ) then
1471  dtype = file_real8
1472  elseif( rp == 4 ) then
1473  dtype = file_real4
1474  else
1475  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
1476  call prc_mpistop
1477  endif
1478  endif
1479 
1480  append_sw = .false.
1481  if ( present(append) ) then
1482  append_sw = append
1483  endif
1484 
1485  call getcftunits(tunits, tsince)
1486 
1487  call filecreate( fid, & ! [OUT]
1488  fileexisted, & ! [OUT]
1489  basename, & ! [IN]
1490  title, & ! [IN]
1491  h_source, & ! [IN]
1492  h_institute, & ! [IN]
1493  prc_masterrank, & ! [IN]
1494  prc_myrank, & ! [IN]
1495  rankidx, & ! [IN]
1496  time_units = tunits, & ! [IN]
1497  append = append_sw ) ! [IN]
1498 
1499  if ( .NOT. fileexisted ) then ! only once
1500  call fileio_set_axes( fid, dtype ) ! [IN]
1501  endif
1502 
1503  if ( axistype == 'ZXYT' ) then
1504  dims = (/'z','x','y'/)
1505  dim1_max = kmax
1506  dim2_max = imaxb
1507  dim3_max = jmaxb
1508  dim1_s = ks
1509  dim1_e = ke
1510  dim2_s = isb
1511  dim2_e = ieb
1512  dim3_s = jsb
1513  dim3_e = jeb
1514  else
1515  write(*,*) 'xxx unsupported axis type. Check!', trim(axistype), ' item:',trim(varname)
1516  call prc_mpistop
1517  endif
1518 
1519  call fileaddvariable( vid, fid, varname, desc, unit, dims, dtype, time_interval ) ! [IN]
1520  allocate( var3d(dim1_max,dim2_max,dim3_max) )
1521 
1522  if ( present(timetarg) ) then
1523  varhalo(:,:,:) = var(:,:,:,timetarg)
1524 
1525  if ( nohalo_ ) then
1526  ! W halo
1527  do k = 1, dim1_max
1528  do j = 1, ja
1529  do i = 1, is-1
1530  varhalo(k,i,j) = rmiss
1531  end do
1532  end do
1533  end do
1534  ! E halo
1535  do k = 1, dim1_max
1536  do j = 1, ja
1537  do i = ie+1, ia
1538  varhalo(k,i,j) = rmiss
1539  end do
1540  end do
1541  end do
1542  ! S halo
1543  do k = 1, dim1_max
1544  do j = 1, js-1
1545  do i = 1, ia
1546  varhalo(k,i,j) = rmiss
1547  end do
1548  end do
1549  end do
1550  ! N halo
1551  do k = 1, dim1_max
1552  do j = je+1, ja
1553  do i = 1, ia
1554  varhalo(k,i,j) = rmiss
1555  end do
1556  end do
1557  end do
1558  end if
1559 
1560  nowtime = (timetarg-1) * time_interval
1561  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
1562  call filewrite( fid, vid, var3d(:,:,:), nowtime, nowtime ) ! [IN]
1563  else
1564  nowtime = 0.0_dp
1565  do n = 1, step
1566  varhalo(:,:,:) = var(:,:,:,n)
1567 
1568  if ( nohalo_ ) then
1569  ! W halo
1570  do k = 1, dim1_max
1571  do j = 1, ja
1572  do i = 1, is-1
1573  varhalo(k,i,j) = rmiss
1574  end do
1575  end do
1576  end do
1577  ! E halo
1578  do k = 1, dim1_max
1579  do j = 1, ja
1580  do i = ie+1, ia
1581  varhalo(k,i,j) = rmiss
1582  end do
1583  end do
1584  end do
1585  ! S halo
1586  do k = 1, dim1_max
1587  do j = 1, js-1
1588  do i = 1, ia
1589  varhalo(k,i,j) = rmiss
1590  end do
1591  end do
1592  end do
1593  ! N halo
1594  do k = 1, dim1_max
1595  do j = je+1, ja
1596  do i = 1, ia
1597  varhalo(k,i,j) = rmiss
1598  end do
1599  end do
1600  end do
1601  end if
1602 
1603  var3d(1:dim1_max,1:dim2_max,1:dim3_max) = varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e)
1604  call filewrite( fid, vid, var3d(:,:,:), nowtime, nowtime ) ! [IN]
1605  nowtime = nowtime + time_interval
1606  enddo
1607  endif
1608 
1609  deallocate( var3d )
1610 
1611  call prof_rapend ('FILE_O_NetCDF', 2)
1612 
1613  return
1614  end subroutine fileio_write_4d
1615 
1616  ! private
1617 
1618  subroutine getcftunits(tunits, date)
1619  character(len=34), intent(out) :: tunits
1620  integer, intent(in) :: date(6)
1621 
1622  write(tunits,'(a,i4.4,"-",i2.2,"-",i2.2," ",i2.2,":",i2.2,":",i2.2)') 'seconds since ', date
1623 
1624  return
1625  end subroutine getcftunits
1626 
1627 end module scale_fileio
1628 !-------------------------------------------------------------------------------
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:150
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:110
real(rp), dimension(:), allocatable, public grid_cxg
center coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fdy
y-length of grid(j+1) to grid(j) [m]
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fxg
face coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fbfy
face buffer factor [0-1]: y
subroutine, public fileio_setup
Setup.
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
module STDIO
Definition: scale_stdio.F90:12
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_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
real(rp), dimension(:), allocatable, public grid_lfz
face coordinate [m]: z, local=global
integer, public ieb
module grid index
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]
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 filesettattr(fid, vname, key, val)
Definition: gtool_file.f90:942
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.
integer, public js
start point of inner domain: y, local
real(rp), dimension(:), allocatable, public grid_cbfx
center buffer factor [0-1]: x
module TIME
Definition: scale_time.F90:15
module PROCESS
real(rp), dimension(:), allocatable, public grid_fbfyg
face buffer factor [0-1]: y, global
real(rp), dimension(:), allocatable, public grid_ucz
center coordinate [m]: z, local=global
module CONSTANT
Definition: scale_const.F90:14
integer, parameter, public prc_masterrank
master process in each communicator
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:), allocatable, public grid_cbfz
center buffer factor [0-1]: z
integer, public prc_myrank
process num in local communicator
real(rp), dimension(:), allocatable, public grid_cx
center coordinate [m]: x, local
module GRID (cartesian)
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
module RM PROCESS
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
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]
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.
subroutine fileio_write_1d(var, basename, title, varname, desc, unit, axistype, datatype, date, subsec, append)
Write 1D data to file.
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