SCALE-RM
Data Types | Functions/Subroutines
scale_fileio Module Reference

module FILE I/O (netcdf) More...

Functions/Subroutines

subroutine, public fileio_setup
 Setup. More...
 
subroutine, public fileio_cleanup
 deallocate buffers More...
 
subroutine, public fileio_set_coordinates (LON, LONX, LONY, LONXY, LAT, LATX, LATY, LATXY, CZ, FZ)
 set latlon and z More...
 
subroutine fileio_check_coordinates_name (basename, atmos, land, urban, transpose)
 check coordinates in the file More...
 
subroutine construct_derived_datatype
 construct MPI derived datatypes for read buffers More...
 
subroutine free_derived_datatype
 free MPI derived datatypes More...
 
subroutine fileio_read_1d (var, basename, varname, axistype, step)
 Read 1D data from file. More...
 
subroutine fileio_read_2d (var, basename, varname, axistype, step)
 Read 2D data from file. More...
 
subroutine fileio_read_3d (var, basename, varname, axistype, step)
 Read 3D data from file. More...
 
subroutine fileio_read_4d (var, basename, varname, axistype, step)
 Read 4D data from file. More...
 
subroutine fileio_read_var_1d (var, fid, varname, axistype, step)
 Read 1D data from file. More...
 
subroutine fileio_read_var_2d (var, fid, varname, axistype, step)
 Read 2D data from file. More...
 
subroutine fileio_read_var_3d (var, fid, varname, axistype, step)
 Read 3D data from file. More...
 
subroutine fileio_read_var_4d (var, fid, varname, axistype, step)
 Read 4D data from file. More...
 
subroutine fileio_write_1d (var, basename, title, varname, desc, unit, axistype, datatype, date, subsec, append)
 Write 1D data to file. More...
 
subroutine fileio_write_2d (var, basename, title, varname, desc, unit, axistype, datatype, date, subsec, append, nohalo, nozcoord)
 Write 2D data to file. More...
 
subroutine fileio_write_3d (var, basename, title, varname, desc, unit, axistype, datatype, date, subsec, append, nohalo)
 Write 3D data to file. More...
 
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. More...
 
subroutine fileio_write_4d (var, basename, title, varname, desc, unit, axistype, datatype, timeintv, tsince, append, timetarg, nohalo)
 Write 4D data to file. More...
 
subroutine, public fileio_open (fid, basename)
 open a netCDF file for read More...
 
subroutine, public fileio_create (fid, basename, title, datatype, date, subsec, append, nozcoord)
 Create/open a netCDF file. More...
 
subroutine, public fileio_enddef (fid)
 Exit netCDF file define mode. More...
 
subroutine, public fileio_flush (fid)
 Flush all pending requests to a netCDF file (PnetCDF only) More...
 
subroutine, public fileio_close (fid)
 Close a netCDF file. More...
 
subroutine fileio_def_axes (fid, dtype, xy)
 define axis variables in the file More...
 
subroutine fileio_write_axes (fid, xy)
 write axis to the file More...
 
subroutine, public fileio_def_var (fid, vid, varname, desc, unit, axistype, datatype, timeintv, nsteps)
 Define a variable to file. More...
 
subroutine fileio_write_var_1d (fid, vid, var, varname, axistype)
 Write 1D data to file. More...
 
subroutine fileio_write_var_2d (fid, vid, var, varname, axistype, nohalo)
 Write 2D data to file. More...
 
subroutine fileio_write_var_3d (fid, vid, var, varname, axistype, nohalo)
 Write 3D data to file. More...
 
subroutine fileio_write_var_3d_t (fid, vid, var, varname, axistype, timeintv, timetarg, nohalo)
 Write 3D data with time dimension to file. More...
 
subroutine fileio_write_var_4d (fid, vid, var, varname, axistype, timeintv, timetarg, nohalo)
 Write 4D data to file. More...
 

Detailed Description

module FILE I/O (netcdf)

Description
general file I/O module frontend interface of netCDF I/O routine
Author
Team SCALE
History
  • 2013-02-20 (H.Yashiro) [new]
NAMELIST
  • PARAM_FILEIO
    nametypedefault valuecomment
    FILEIO_DATACHECK_CRITERIA real(RP)

History Output
No history output

Function/Subroutine Documentation

◆ fileio_setup()

subroutine, public scale_fileio::fileio_setup ( )

Setup.

Definition at line 145 of file scale_fileio.F90.

References scale_const::const_eps, construct_derived_datatype(), scale_stdio::h_institute, scale_stdio::h_source, scale_grid_index::ia, scale_grid_index::ieb, scale_grid_index::ihalo, scale_grid_index::imax, scale_stdio::io_aggregate, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_grid_index::isb, scale_grid_index::ja, scale_grid_index::jeb, scale_grid_index::jhalo, scale_grid_index::jmax, scale_grid_index::jsb, scale_grid_index::ka, scale_rm_process::prc_2drank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_rm_process::prc_num_x, and scale_rm_process::prc_num_y.

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

145  use scale_process, only: &
146  prc_myrank, &
148  use scale_rm_process, only: &
149  prc_2drank, &
150  prc_num_x, &
151  prc_num_y
152  use scale_const, only: &
153  eps => const_eps
154  implicit none
155 
156  namelist / param_fileio / &
157  fileio_datacheck_criteria
158 
159  integer :: ierr
160  !---------------------------------------------------------------------------
161 
162  if( io_l ) write(io_fid_log,*)
163  if( io_l ) write(io_fid_log,*) '++++++ Module[FIELIO] / Categ[ATMOS-RM IO] / Origin[SCALElib]'
164 
165  fileio_datacheck_criteria = eps * 10.0_rp
166 
167  !--- read namelist
168  rewind(io_fid_conf)
169  read(io_fid_conf,nml=param_fileio,iostat=ierr)
170  if( ierr < 0 ) then !--- missing
171  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
172  elseif( ierr > 0 ) then !--- fatal error
173  write(*,*) 'xxx Not appropriate names in namelist PARAM_FILEIO. Check!'
174  call prc_mpistop
175  endif
176  if( io_nml ) write(io_fid_nml,nml=param_fileio)
177 
178  if( io_l ) write(io_fid_log,*)
179  if( io_l ) write(io_fid_log,*) '*** NetCDF header information ***'
180  if( io_l ) write(io_fid_log,*) '*** Data source : ', trim(h_source)
181  if( io_l ) write(io_fid_log,*) '*** Institute : ', trim(h_institute)
182  if( io_l ) write(io_fid_log,*)
183  if( io_l ) write(io_fid_log,*) '*** Data consistency criteria : ', &
184  '(file-internal)/internal = ', fileio_datacheck_criteria
185 
186  if ( io_aggregate ) then
187  ! construct indices independent from PRC_PERIODIC_X/Y
188  xsb = 1 + ihalo
189  if( prc_2drank(prc_myrank,1) == 0 ) xsb = 1
190  xeb = imax + ihalo
191  if( prc_2drank(prc_myrank,1) == prc_num_x-1 ) xeb = ia
192 
193  ysb = 1 + jhalo
194  if( prc_2drank(prc_myrank,2) == 0 ) ysb = 1
195  yeb = jmax + jhalo
196  if( prc_2drank(prc_myrank,2) == prc_num_y-1 ) yeb = ja
197  else
198  xsb = isb
199  xeb = ieb
200  ysb = jsb
201  yeb = jeb
202  end if
203 
204  allocate( axis_lon( ia, ja) )
205  allocate( axis_lonx(0:ia, ja) )
206  allocate( axis_lony( ia,0:ja) )
207  allocate( axis_lonxy(0:ia,0:ja) )
208  allocate( axis_lat( ia, ja) )
209  allocate( axis_latx(0:ia, ja) )
210  allocate( axis_laty( ia,0:ja) )
211  allocate( axis_latxy(0:ia,0:ja) )
212 
213  allocate( axis_hzxy( ka,ia,ja) )
214  allocate( axis_hwxy(0:ka,ia,ja) )
215 
216  if( io_aggregate ) call construct_derived_datatype
217 
218  file_closed(:) = .true.
219 
220  return
integer, public imax
of computational cells: x, local
integer, public prc_num_x
x length of 2D processor topology
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
integer, public jhalo
of halo cells: y
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public prc_myrank
process num in local communicator
module RM PROCESS
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, public isb
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public jsb
integer, public jmax
of computational cells: y, local
integer, public ihalo
of halo cells: x
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_cleanup()

subroutine, public scale_fileio::fileio_cleanup ( )

deallocate buffers

Definition at line 226 of file scale_fileio.F90.

References free_derived_datatype().

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

226  implicit none
227  !---------------------------------------------------------------------------
228 
229  deallocate( axis_lon )
230  deallocate( axis_lonx )
231  deallocate( axis_lony )
232  deallocate( axis_lonxy )
233  deallocate( axis_lat )
234  deallocate( axis_latx )
235  deallocate( axis_laty )
236  deallocate( axis_latxy )
237  deallocate( axis_hzxy )
238  deallocate( axis_hwxy )
239 
240  call free_derived_datatype
241 
242  call closeall
243 
244  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_set_coordinates()

subroutine, public scale_fileio::fileio_set_coordinates ( real(rp), dimension ( ia, ja), intent(in)  LON,
real(rp), dimension (0:ia, ja), intent(in)  LONX,
real(rp), dimension ( ia,0:ja), intent(in)  LONY,
real(rp), dimension(0:ia,0:ja), intent(in)  LONXY,
real(rp), dimension ( ia, ja), intent(in)  LAT,
real(rp), dimension (0:ia, ja), intent(in)  LATX,
real(rp), dimension ( ia,0:ja), intent(in)  LATY,
real(rp), dimension(0:ia,0:ja), intent(in)  LATXY,
real(rp), dimension ( ka,ia,ja), intent(in)  CZ,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ 
)

set latlon and z

Definition at line 260 of file scale_fileio.F90.

References scale_const::const_d2r.

Referenced by scale_grid_real::real_setup(), and scale_grid_real::real_update_z().

260  use scale_const, only: &
261  d2r => const_d2r
262  implicit none
263 
264  real(RP), intent(in) :: LON ( IA, JA)
265  real(RP), intent(in) :: LONX (0:IA, JA)
266  real(RP), intent(in) :: LONY ( IA,0:JA)
267  real(RP), intent(in) :: LONXY(0:IA,0:JA)
268  real(RP), intent(in) :: LAT ( IA, JA)
269  real(RP), intent(in) :: LATX (0:IA, JA)
270  real(RP), intent(in) :: LATY ( IA,0:JA)
271  real(RP), intent(in) :: LATXY(0:IA,0:JA)
272  real(RP), intent(in) :: CZ ( KA,IA,JA)
273  real(RP), intent(in) :: FZ (0:KA,IA,JA)
274  !---------------------------------------------------------------------------
275 
276  axis_lon(:,:) = lon(:,:) / d2r
277  axis_lonx(:,:) = lonx(:,:) / d2r
278  axis_lony(:,:) = lony(:,:) / d2r
279  axis_lonxy(:,:) = lonxy(:,:) / d2r
280  axis_lat(:,:) = lat(:,:) / d2r
281  axis_latx(:,:) = latx(:,:) / d2r
282  axis_laty(:,:) = laty(:,:) / d2r
283  axis_latxy(:,:) = latxy(:,:) / d2r
284 
285  axis_hzxy(:,:,:) = cz(:,:,:)
286  axis_hwxy(:,:,:) = fz(:,:,:)
287 
288  set_coordinates = .true.
289 
290  return
module CONSTANT
Definition: scale_const.F90:14
Here is the caller graph for this function:

◆ fileio_check_coordinates_name()

subroutine scale_fileio::fileio_check_coordinates_name ( character(len=*), intent(in)  basename,
logical, intent(in), optional  atmos,
logical, intent(in), optional  land,
logical, intent(in), optional  urban,
logical, intent(in), optional  transpose 
)

check coordinates in the file

Parameters
[in]basenamebasename of the file
[in]atmoscheck atmospheric coordinates
[in]landcheck land coordinates
[in]urbancheck urban coordinates

Definition at line 301 of file scale_fileio.F90.

References fileio_flush(), fileio_open(), fileio_read_var_1d(), fileio_read_var_2d(), fileio_read_var_3d(), scale_grid::grid_cx, scale_grid::grid_cy, scale_grid::grid_cz, scale_land_grid::grid_lcz, scale_urban_grid::grid_ucz, scale_grid_index::ieb, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, scale_land_grid_index::lke, scale_land_grid_index::lks, scale_urban_grid_index::uke, and scale_urban_grid_index::uks.

301  implicit none
302 
303  character(len=*), intent(in) :: basename
304  logical, intent(in), optional :: atmos
305  logical, intent(in), optional :: land
306  logical, intent(in), optional :: urban
307  logical, intent(in), optional :: transpose
308 
309  logical :: atmos_
310  logical :: land_
311  logical :: urban_
312  logical :: transpose_
313 
314  integer :: fid
315  !---------------------------------------------------------------------------
316 
317  atmos_ = .false.
318  land_ = .false.
319  urban_ = .false.
320  transpose_ = .false.
321 
322  if( present(atmos) ) atmos_ = atmos
323  if( present(land ) ) land_ = land
324  if( present(urban) ) urban_ = urban
325  if( present(transpose) ) transpose_ = transpose
326 
327  call fileio_open( fid, & ! [OUT]
328  basename ) ! [IN]
329 
330  call fileio_check_coordinates_id( fid, & ! [IN]
331  atmos_, land_, urban_, & ! [IN]
332  transpose_ )
333 
334  return
Here is the call graph for this function:

◆ construct_derived_datatype()

subroutine scale_fileio::construct_derived_datatype ( )

construct MPI derived datatypes for read buffers

Definition at line 434 of file scale_fileio.F90.

References scale_grid_index::ia, scale_grid_index::ihalo, scale_grid_index::imaxb, scale_grid_index::is_ing, scale_grid_index::isb, scale_grid_index::ja, scale_grid_index::jhalo, scale_grid_index::js_ing, scale_grid_index::ka, scale_grid_index::khalo, scale_grid_index::kmax, scale_grid_index::ks, scale_land_grid_index::lkmax, scale_land_grid_index::lks, scale_process::prc_mpistop(), scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_precision::rp, scale_urban_grid_index::ukmax, and scale_urban_grid_index::uks.

Referenced by fileio_setup().

434  use mpi
435  use scale_process, only: &
437  use scale_rm_process, only: &
438  prc_num_x, &
439  prc_num_y
440  implicit none
441 
442  integer :: err, order
443  integer :: sizes(3), subsizes(3), sub_off(3)
444  !---------------------------------------------------------------------------
445 
446  order = mpi_order_fortran
447 
448  centertypexy = mpi_datatype_null
449  centertypezx = mpi_datatype_null
450  centertypezxy = mpi_datatype_null
451  centertypeland = mpi_datatype_null
452  centertypeurban = mpi_datatype_null
453 
454  etype = mpi_float
455 
456  if( rp == 8 ) etype = mpi_double_precision
457 
458  ! for axistype == 'XY'
459  startxy(1) = is_ing - ihalo
460  startxy(2) = js_ing - jhalo
461  countxy(1) = ia
462  countxy(2) = ja
463  ! for axistype == 'ZXY'
464  startzxy(1) = 1
465  startzxy(2:3) = startxy(1:2)
466  countzxy(1) = kmax
467  countzxy(2:3) = countxy(1:2)
468  ! construct MPI subarray data type
469  sizes(1) = ka
470  sizes(2) = ia
471  sizes(3) = ja
472  subsizes(1) = kmax
473  subsizes(2) = ia
474  subsizes(3) = ja
475  sub_off(1) = ks - 1 ! MPI start index starts with 0
476  sub_off(2) = 0
477  sub_off(3) = 0
478 
479  call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypezxy, err)
480  call mpi_type_commit(centertypezxy, err)
481 
482  ! for axistype == 'Land'
483  startland(1) = 1
484  startland(2:3) = startxy(1:2)
485  countland(1) = lkmax
486  countland(2:3) = countxy(1:2)
487  ! construct MPI subarray data type
488  sizes(1) = lkmax
489  subsizes(1) = lkmax
490  sub_off(1) = lks - 1 ! MPI start index starts with 0
491 
492  call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypeland, err)
493  call mpi_type_commit(centertypeland, err)
494 
495  ! for axistype == 'URBAN'
496  starturban(1) = 1
497  starturban(2:3) = startxy(1:2)
498  counturban(1) = ukmax
499  counturban(2:3) = countxy(1:2)
500  ! construct MPI subarray data type
501  sizes(1) = ukmax
502  subsizes(1) = ukmax
503  sub_off(1) = uks - 1 ! MPI start index starts with 0
504 
505  call mpi_type_create_subarray(3, sizes, subsizes, sub_off, order, etype, centertypeurban, err)
506  call mpi_type_commit(centertypeurban, err)
507 
508  ! for axistype == 'ZX'
509  startzx(1) = khalo+1
510  startzx(2) = is_ing - ihalo
511  countzx(1) = khalo
512  countzx(2) = ia
513  ! construct MPI subarray data type
514  sizes(1) = ka
515  sizes(2) = ia
516  subsizes(1) = kmax
517  subsizes(2) = imaxb
518  sub_off(1) = khalo ! MPI start index starts with 0
519  sub_off(2) = isb - 1 ! MPI start index starts with 0
520 
521  call mpi_type_create_subarray(2, sizes, subsizes, sub_off, order, etype, centertypezx, err)
522  call mpi_type_commit(centertypezx, err)
523 
524  return
integer, public prc_num_x
x length of 2D processor topology
subroutine, public prc_mpistop
Abort MPI.
integer, public imaxb
integer, public prc_num_y
y length of 2D processor topology
integer, public jhalo
of halo cells: y
module PROCESS
integer, public is_ing
start point of the inner domain: cx, global
integer, parameter, public khalo
of halo cells: z
module RM PROCESS
integer, public js_ing
start point of the inner domain: cy, global
integer, public isb
integer, public ihalo
of halo cells: x
Here is the call graph for this function:
Here is the caller graph for this function:

◆ free_derived_datatype()

subroutine scale_fileio::free_derived_datatype ( )

free MPI derived datatypes

Definition at line 530 of file scale_fileio.F90.

References scale_stdio::io_aggregate.

Referenced by fileio_cleanup().

530  use mpi
531  implicit none
532 
533  integer :: err
534  !---------------------------------------------------------------------------
535 
536  if( .NOT. io_aggregate ) return
537 
538  if( centertypexy /= mpi_datatype_null ) call mpi_type_free(centertypexy, err)
539  if( centertypezx /= mpi_datatype_null ) call mpi_type_free(centertypezx, err)
540  if( centertypezxy /= mpi_datatype_null ) call mpi_type_free(centertypezxy, err)
541  if( centertypeland /= mpi_datatype_null ) call mpi_type_free(centertypeland, err)
542  if( centertypeurban /= mpi_datatype_null ) call mpi_type_free(centertypeurban, err)
543 
544  return
Here is the caller graph for this function:

◆ fileio_read_1d()

subroutine scale_fileio::fileio_read_1d ( real(rp), dimension(:), intent(out)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 1D data from file.

Parameters
[out]varvalue of the variable
[in]basenamebasename of the file
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y)
[in]stepstep number

Definition at line 555 of file scale_fileio.F90.

References fileio_close(), fileio_open(), fileio_read_var_1d(), and scale_process::prc_mpistop().

555  use gtool_file, only: &
556  fileread
557  use scale_process, only: &
559  implicit none
560 
561  real(RP), intent(out) :: var(:)
562  character(len=*), intent(in) :: basename
563  character(len=*), intent(in) :: varname
564  character(len=*), intent(in) :: axistype
565  integer, intent(in) :: step
566 
567  integer :: fid
568  !---------------------------------------------------------------------------
569 
570  call fileio_open( fid, & ! [OUT]
571  basename ) ! [IN]
572 
573  call fileio_read_var_1d( var(:), & ! [OUT]
574  fid, varname, axistype, step ) ! [IN]
575 
576  call fileio_close( fid )
577 
578  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:

◆ fileio_read_2d()

subroutine scale_fileio::fileio_read_2d ( real(rp), dimension(:,:), intent(out)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 2D data from file.

Parameters
[out]varvalue of the variable
[in]basenamebasename of the file
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y)
[in]stepstep number

Definition at line 589 of file scale_fileio.F90.

References fileio_close(), fileio_open(), fileio_read_var_2d(), and scale_process::prc_mpistop().

589  use gtool_file, only: &
590  fileread
591  use scale_process, only: &
593  implicit none
594 
595  real(RP), intent(out) :: var(:,:)
596  character(len=*), intent(in) :: basename
597  character(len=*), intent(in) :: varname
598  character(len=*), intent(in) :: axistype
599  integer, intent(in) :: step
600 
601  integer :: fid
602  !---------------------------------------------------------------------------
603 
604  call fileio_open( fid, & ! [OUT]
605  basename ) ! [IN]
606 
607  call fileio_read_var_2d( var(:,:), & ! [OUT]
608  fid, varname, axistype, step ) ! [IN]
609 
610  call fileio_close( fid )
611 
612  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:

◆ fileio_read_3d()

subroutine scale_fileio::fileio_read_3d ( real(rp), dimension(:,:,:), intent(out)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 3D data from file.

Parameters
[out]varvalue of the variable
[in]basenamebasename of the file
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y/T)
[in]stepstep number

Definition at line 623 of file scale_fileio.F90.

References fileio_close(), fileio_open(), fileio_read_var_3d(), and scale_process::prc_mpistop().

623  use gtool_file, only: &
624  fileread
625  use scale_process, only: &
627  implicit none
628 
629  real(RP), intent(out) :: var(:,:,:)
630  character(len=*), intent(in) :: basename
631  character(len=*), intent(in) :: varname
632  character(len=*), intent(in) :: axistype
633  integer, intent(in) :: step
634 
635  integer :: fid
636  !---------------------------------------------------------------------------
637 
638  call fileio_open( fid, & ! [OUT]
639  basename ) ! [IN]
640 
641  call fileio_read_var_3d( var(:,:,:), & ! [OUT]
642  fid, varname, axistype, step ) ! [IN]
643 
644  call fileio_close( fid )
645 
646  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:

◆ fileio_read_4d()

subroutine scale_fileio::fileio_read_4d ( real(rp), dimension(:,:,:,:), intent(out)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 4D data from file.

Parameters
[out]varvalue of the variable
[in]basenamebasename of the file
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y/Time)
[in]stepstep number

Definition at line 657 of file scale_fileio.F90.

References fileio_close(), fileio_open(), fileio_read_var_4d(), and scale_process::prc_mpistop().

657  use gtool_file, only: &
658  fileread
659  use scale_process, only: &
661  implicit none
662 
663  real(RP), intent(out) :: var(:,:,:,:)
664  character(len=*), intent(in) :: basename
665  character(len=*), intent(in) :: varname
666  character(len=*), intent(in) :: axistype
667  integer, intent(in) :: step
668 
669  integer :: fid
670  !---------------------------------------------------------------------------
671 
672  call fileio_open( fid, & ! [OUT]
673  basename ) ! [IN]
674 
675  call fileio_read_var_4d( var(:,:,:,:), & ! [OUT]
676  fid, varname, axistype, step ) ! [IN]
677 
678  call fileio_close( fid )
679 
680  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:

◆ fileio_read_var_1d()

subroutine scale_fileio::fileio_read_var_1d ( real(rp), dimension(:), intent(out)  var,
integer, intent(in)  fid,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 1D data from file.

Parameters
[out]varvalue of the variable
[in]fidfile ID
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y)
[in]stepstep number

Definition at line 691 of file scale_fileio.F90.

References scale_grid_index::ia, scale_grid_index::ieb, scale_grid_index::ihalo, scale_stdio::io_aggregate, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is_ing, scale_grid_index::isb, scale_grid_index::ja, scale_grid_index::jeb, scale_grid_index::jhalo, scale_grid_index::js_ing, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, scale_land_grid_index::lkmax, scale_process::prc_mpistop(), scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_urban_grid_index::ukmax.

Referenced by fileio_check_coordinates_name(), and fileio_read_1d().

691  use gtool_file, only: &
692  fileread
693  use scale_process, only: &
695  use scale_rm_process, only: &
696  prc_num_x, &
697  prc_num_y
698  use mpi
699  implicit none
700 
701  real(RP), intent(out) :: var(:)
702  integer, intent(in) :: fid
703  character(len=*), intent(in) :: varname
704  character(len=*), intent(in) :: axistype
705  integer, intent(in) :: step
706 
707  integer :: dim1_S, dim1_E
708  integer :: start(1) ! start offset of globale variable
709  integer :: count(1) ! request length to the global variable
710  !---------------------------------------------------------------------------
711 
712  call prof_rapstart('FILE_I_NetCDF', 2)
713 
714  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Read from file (1D), name : ', trim(varname)
715 
716  if ( io_aggregate ) then
717  ! read data and halos into the local buffer
718  if ( axistype == 'Z' ) then
719  start(1) = 1
720  count(1) = kmax
721  call fileread( var(ks:ke), fid, varname, step, &
722  ntypes=kmax, dtype=etype, start=start, count=count )
723  elseif( axistype == 'LZ' ) then
724  start(1) = 1
725  count(1) = lkmax
726  call fileread( var, fid, varname, step, &
727  ntypes=lkmax, dtype=etype, start=start, count=count )
728  elseif( axistype == 'UZ' ) then
729  start(1) = 1
730  count(1) = ukmax
731  call fileread( var, fid, varname, step, &
732  ntypes=ukmax, dtype=etype, start=start, count=count )
733  elseif( axistype == 'X' .OR. axistype == 'CX' ) then
734  start(1) = is_ing - ihalo
735  count(1) = ia
736  call fileread( var, fid, varname, step, &
737  ntypes=ia, dtype=etype, start=start, count=count )
738  elseif( axistype == 'Y' .OR. axistype == 'CY' ) then
739  start(1) = js_ing - jhalo
740  count(1) = ja
741  call fileread( var, fid, varname, step, &
742  ntypes=ja, dtype=etype, start=start, count=count )
743  else
744  write(*,*) 'xxx [FILEIO_read_var_1D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
745  call prc_mpistop
746  endif
747  else
748  if ( axistype == 'Z' ) then
749  dim1_s = ks
750  dim1_e = ke
751  elseif( axistype == 'LZ' ) then
752  dim1_s = 1
753  dim1_e = lkmax
754  elseif( axistype == 'UZ' ) then
755  dim1_s = 1
756  dim1_e = ukmax
757  elseif( axistype == 'X' ) then
758  dim1_s = isb
759  dim1_e = ieb
760  elseif( axistype == 'CX' ) then
761  dim1_s = 1
762  dim1_e = ia
763  elseif( axistype == 'Y' ) then
764  dim1_s = jsb
765  dim1_e = jeb
766  elseif( axistype == 'CY' ) then
767  dim1_s = 1
768  dim1_e = ja
769  else
770  write(*,*) 'xxx [FILEIO_read_var_1D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
771  call prc_mpistop
772  endif
773 
774  call fileread( var(dim1_s:dim1_e), fid, varname, step )
775  end if
776 
777  call prof_rapend ('FILE_I_NetCDF', 2)
778 
779  return
integer, public prc_num_x
x length of 2D processor topology
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
integer, public jhalo
of halo cells: y
module PROCESS
integer, public is_ing
start point of the inner domain: cx, global
module RM PROCESS
integer, public js_ing
start point of the inner domain: cy, global
integer, public isb
integer, public jsb
integer, public ihalo
of halo cells: x
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_read_var_2d()

subroutine scale_fileio::fileio_read_var_2d ( real(rp), dimension(:,:), intent(out)  var,
integer, intent(in)  fid,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 2D data from file.

Parameters
[out]varvalue of the variable
[in]fidfile ID
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y)
[in]stepstep number

Definition at line 790 of file scale_fileio.F90.

References scale_grid_index::ia, scale_grid_index::ieb, scale_stdio::io_aggregate, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::isb, scale_grid_index::ja, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, scale_process::prc_mpistop(), scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by fileio_check_coordinates_name(), and fileio_read_2d().

790  use gtool_file, only: &
791  fileread
792  use scale_process, only: &
794  use scale_rm_process, only: &
795  prc_num_x, &
796  prc_num_y
797  use mpi
798  implicit none
799 
800  real(RP), intent(out) :: var(:,:)
801  integer, intent(in) :: fid
802  character(len=*), intent(in) :: varname
803  character(len=*), intent(in) :: axistype
804  integer, intent(in) :: step
805 
806  integer :: dim1_S, dim1_E
807  integer :: dim2_S, dim2_E
808  !---------------------------------------------------------------------------
809 
810  call prof_rapstart('FILE_I_NetCDF', 2)
811 
812  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Read from file (2D), name : ', trim(varname)
813 
814  if ( io_aggregate ) then
815  ! read data and halos into the local buffer
816  if ( axistype == 'XY' ) then
817  call fileread( var, fid, varname, step, &
818  ntypes=ia*ja, dtype=etype, start=startxy, count=countxy )
819  elseif( axistype == 'ZX' ) then
820  ! Because KHALO is not saved in files, we use centerTypeZX, an MPI
821  ! derived datatype to describe the layout of local read buffer
822  call fileread( var, fid, varname, step, &
823  ntypes=1, dtype=centertypezx, start=startzx, count=countzx )
824  else
825  write(*,*) 'xxx [FILEIO_read_var_2D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
826  call prc_mpistop
827  endif
828  else
829  if ( axistype == 'XY' ) then
830  dim1_s = isb
831  dim1_e = ieb
832  dim2_s = jsb
833  dim2_e = jeb
834  elseif( axistype == 'ZX' ) then
835  dim1_s = ks
836  dim1_e = ke
837  dim2_s = isb
838  dim2_e = ieb
839  else
840  write(*,*) 'xxx [FILEIO_read_var_2D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
841  call prc_mpistop
842  endif
843 
844  call fileread( var(dim1_s:dim1_e,dim2_s:dim2_e), fid, varname, step )
845  end if
846 
847  call prof_rapend ('FILE_I_NetCDF', 2)
848 
849  return
integer, public prc_num_x
x length of 2D processor topology
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
module PROCESS
module RM PROCESS
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_read_var_3d()

subroutine scale_fileio::fileio_read_var_3d ( real(rp), dimension(:,:,:), intent(out)  var,
integer, intent(in)  fid,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 3D data from file.

Parameters
[out]varvalue of the variable
[in]fidfile ID
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y/T)
[in]stepstep number

Definition at line 860 of file scale_fileio.F90.

References scale_grid_index::ia, scale_grid_index::ieb, scale_stdio::io_aggregate, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::isb, scale_grid_index::ja, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, scale_land_grid_index::lke, scale_land_grid_index::lks, scale_process::prc_mpistop(), scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_urban_grid_index::uke, and scale_urban_grid_index::uks.

Referenced by fileio_check_coordinates_name(), and fileio_read_3d().

860  use gtool_file, only: &
861  fileread
862  use scale_process, only: &
864  use scale_rm_process, only: &
865  prc_num_x, &
866  prc_num_y
867  implicit none
868 
869  real(RP), intent(out) :: var(:,:,:)
870  integer, intent(in) :: fid
871  character(len=*), intent(in) :: varname
872  character(len=*), intent(in) :: axistype
873  integer, intent(in) :: step
874 
875  integer :: dim1_S, dim1_E
876  integer :: dim2_S, dim2_E
877  integer :: dim3_S, dim3_E
878  !---------------------------------------------------------------------------
879 
880  call prof_rapstart('FILE_I_NetCDF', 2)
881 
882  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Read from file (3D), name : ', trim(varname)
883 
884  if ( io_aggregate ) then
885  ! read data and halos into the local buffer
886  ! Because KHALO is not saved in files, we use mpi derived datatypes to
887  ! describe the layout of local read buffer
888  if ( axistype == 'ZXY' ) then
889  call fileread( var, fid, varname, step, &
890  ntypes=1, dtype=centertypezxy, start=startzxy, count=countzxy )
891  elseif( axistype == 'XYT' ) then
892  startxy(3) = 1
893  countxy(3) = step
894  call fileread( var, fid, varname, step, &
895  ntypes=step*ia*ja, dtype=etype, start=startxy, count=countxy )
896  elseif( axistype == 'Land' ) then
897  call fileread( var, fid, varname, step, &
898  ntypes=1, dtype=centertypeland, start=startland, count=countland )
899  elseif( axistype == 'Urban' ) then
900  call fileread( var, fid, varname, step, &
901  ntypes=1, dtype=centertypeurban, start=starturban, count=counturban )
902  else
903  write(*,*) 'xxx [FILEIO_read_var_3D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
904  call prc_mpistop
905  endif
906  else
907  if ( axistype == 'ZXY' ) then
908  dim1_s = ks
909  dim1_e = ke
910  dim2_s = isb
911  dim2_e = ieb
912  dim3_s = jsb
913  dim3_e = jeb
914  elseif( axistype == 'XYT' ) then
915  dim1_s = isb
916  dim1_e = ieb
917  dim2_s = jsb
918  dim2_e = jeb
919  dim3_s = 1
920  dim3_e = step
921  elseif( axistype == 'Land' ) then
922  dim1_s = lks
923  dim1_e = lke
924  dim2_s = isb
925  dim2_e = ieb
926  dim3_s = jsb
927  dim3_e = jeb
928  elseif( axistype == 'Urban' ) then
929  dim1_s = uks
930  dim1_e = uke
931  dim2_s = isb
932  dim2_e = ieb
933  dim3_s = jsb
934  dim3_e = jeb
935  else
936  write(*,*) 'xxx [FILEIO_read_var_3D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
937  call prc_mpistop
938  endif
939 
940  call fileread( var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e), &
941  fid, varname, step )
942  end if
943 
944  call prof_rapend ('FILE_I_NetCDF', 2)
945 
946  return
integer, public prc_num_x
x length of 2D processor topology
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
module PROCESS
module RM PROCESS
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_read_var_4d()

subroutine scale_fileio::fileio_read_var_4d ( real(rp), dimension(:,:,:,:), intent(out)  var,
integer, intent(in)  fid,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
integer, intent(in)  step 
)

Read 4D data from file.

Parameters
[out]varvalue of the variable
[in]fidfile ID
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y/Time)
[in]stepstep number

Definition at line 957 of file scale_fileio.F90.

References scale_grid_index::ieb, scale_stdio::io_aggregate, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, scale_process::prc_mpistop(), scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by fileio_read_4d().

957  use gtool_file, only: &
958  fileread
959  use scale_process, only: &
961  use scale_rm_process, only: &
962  prc_num_x, &
963  prc_num_y
964  implicit none
965 
966  real(RP), intent(out) :: var(:,:,:,:)
967  integer, intent(in) :: fid
968  character(len=*), intent(in) :: varname
969  character(len=*), intent(in) :: axistype
970  integer, intent(in) :: step
971 
972  integer :: dim1_S, dim1_E
973  integer :: dim2_S, dim2_E
974  integer :: dim3_S, dim3_E
975  integer :: dim4_S, dim4_E
976  !---------------------------------------------------------------------------
977 
978  call prof_rapstart('FILE_I_NetCDF', 2)
979 
980  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Read from file (4D), name : ', trim(varname)
981 
982  if ( io_aggregate ) then
983  ! read data and halos into the local buffer
984  if ( axistype == 'ZXYT' ) then
985  startzxy(4) = 1
986  countzxy(4) = step
987  call fileread( var, fid, varname, step, &
988  ntypes=step, dtype=centertypezxy, start=startzxy, count=countzxy )
989  else
990  write(*,*) 'xxx [FILEIO_read_var_4D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
991  call prc_mpistop
992  endif
993  else
994  if ( axistype == 'ZXYT' ) then
995  dim1_s = ks
996  dim1_e = ke
997  dim2_s = isb
998  dim2_e = ieb
999  dim3_s = jsb
1000  dim3_e = jeb
1001  dim4_s = 1
1002  dim4_e = step
1003  else
1004  write(*,*) 'xxx [FILEIO_read_var_4D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
1005  call prc_mpistop
1006  endif
1007 
1008  call fileread( var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e,dim4_s:dim4_e), &
1009  fid, varname, step )
1010  end if
1011 
1012  call prof_rapend ('FILE_I_NetCDF', 2)
1013 
1014  return
integer, public prc_num_x
x length of 2D processor topology
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
module PROCESS
module RM PROCESS
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_write_1d()

subroutine scale_fileio::fileio_write_1d ( real(rp), dimension(:), intent(in)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  title,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  desc,
character(len=*), intent(in)  unit,
character(len=*), intent(in)  axistype,
character(len=*), intent(in)  datatype,
integer, dimension(6), intent(in), optional  date,
real(dp), intent(in), optional  subsec,
logical, intent(in), optional  append 
)

Write 1D data to file.

Parameters
[in]varvalue of the variable
[in]basenamebasename of the file
[in]titletitle of the file
[in]varnamename of the variable
[in]descdescription of the variable
[in]unitunit of the variable
[in]axistypeaxis type (Z/X/Y)
[in]datatypedata type (REAL8/REAL4/default)
[in]dateymdhms of the time
[in]subsecsubsec of the time
[in]appendswitch whether append existing file or not (default=false)

Definition at line 1031 of file scale_fileio.F90.

References gtool_file::filecreate(), fileio_create(), fileio_def_var(), fileio_enddef(), fileio_write_var_1d(), scale_stdio::io_fid_log, scale_stdio::io_l, and scale_process::prc_mpistop().

1031  use gtool_file, only: &
1032  filecreate, &
1033  fileaddvariable, &
1034  filesetglobalattribute, &
1035  filewrite
1036  use scale_process, only: &
1037  prc_mpistop
1038  implicit none
1039 
1040  real(RP), intent(in) :: var(:)
1041  character(len=*), intent(in) :: basename
1042  character(len=*), intent(in) :: title
1043  character(len=*), intent(in) :: varname
1044  character(len=*), intent(in) :: desc
1045  character(len=*), intent(in) :: unit
1046  character(len=*), intent(in) :: axistype
1047  character(len=*), intent(in) :: datatype
1048 
1049  integer, intent(in), optional :: date(6)
1050  real(DP), intent(in), optional :: subsec
1051  logical, intent(in), optional :: append
1052 
1053  integer :: fid, vid
1054  !---------------------------------------------------------------------------
1055 
1056  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Write to file (1D), name : ', trim(varname)
1057 
1058  call fileio_create( fid, & ! [OUT]
1059  basename, title, datatype, date, subsec, append )
1060 
1061  call fileio_def_var( fid, vid, varname, desc, unit, axistype, datatype )
1062 
1063  call fileio_enddef( fid )
1064 
1065  call fileio_write_var_1d( fid, vid, var, varname, axistype )
1066 
1067  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
subroutine, public filecreate(fid, existed, basename, title, source, institution, master, myrank, rankidx, single, time_units, append, comm)
Definition: gtool_file.f90:191
Here is the call graph for this function:

◆ fileio_write_2d()

subroutine scale_fileio::fileio_write_2d ( real(rp), dimension(:,:), intent(in)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  title,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  desc,
character(len=*), intent(in)  unit,
character(len=*), intent(in)  axistype,
character(len=*), intent(in)  datatype,
integer, dimension(6), intent(in), optional  date,
real(dp), intent(in), optional  subsec,
logical, intent(in), optional  append,
logical, intent(in), optional  nohalo,
logical, intent(in), optional  nozcoord 
)

Write 2D data to file.

Parameters
[in]varvalue of the variable
[in]basenamebasename of the file
[in]titletitle of the file
[in]varnamename of the variable
[in]descdescription of the variable
[in]unitunit of the variable
[in]axistypeaxis type (Z/X/Y)
[in]datatypedata type (REAL8/REAL4/default)
[in]dateymdhms of the time
[in]subsecsubsec of the time
[in]appendswitch whether append existing file or not (default=false)
[in]nohaloswitch whether include halo data or not (default=false)
[in]nozcoordswitch whether include zcoordinate or not (default=false)

Definition at line 1086 of file scale_fileio.F90.

References gtool_file::filecreate(), fileio_create(), fileio_def_var(), fileio_enddef(), fileio_write_var_2d(), scale_stdio::io_fid_log, scale_stdio::io_l, and scale_process::prc_mpistop().

1086  use gtool_file, only: &
1087  filecreate, &
1088  fileaddvariable, &
1089  filesetglobalattribute, &
1090  filewrite
1091  use scale_process, only: &
1092  prc_mpistop
1093  implicit none
1094 
1095  real(RP), intent(in) :: var(:,:)
1096  character(len=*), intent(in) :: basename
1097  character(len=*), intent(in) :: title
1098  character(len=*), intent(in) :: varname
1099  character(len=*), intent(in) :: desc
1100  character(len=*), intent(in) :: unit
1101  character(len=*), intent(in) :: axistype
1102  character(len=*), intent(in) :: datatype
1103 
1104  integer, intent(in), optional :: date(6)
1105  real(DP), intent(in), optional :: subsec
1106  logical, intent(in), optional :: append
1107  logical, intent(in), optional :: nohalo
1108  logical, intent(in), optional :: nozcoord
1109 
1110  integer :: fid, vid
1111  !---------------------------------------------------------------------------
1112 
1113  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Write to file (2D), name : ', trim(varname)
1114 
1115  call fileio_create( fid, & ! [OUT]
1116  basename, title, datatype, date, subsec, append, nozcoord )
1117 
1118  call fileio_def_var( fid, vid, varname, desc, unit, axistype, datatype )
1119 
1120  call fileio_enddef( fid )
1121 
1122  call fileio_write_var_2d( fid, vid, var, varname, axistype, nohalo )
1123 
1124  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
subroutine, public filecreate(fid, existed, basename, title, source, institution, master, myrank, rankidx, single, time_units, append, comm)
Definition: gtool_file.f90:191
Here is the call graph for this function:

◆ fileio_write_3d()

subroutine scale_fileio::fileio_write_3d ( real(rp), dimension(:,:,:), intent(in)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  title,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  desc,
character(len=*), intent(in)  unit,
character(len=*), intent(in)  axistype,
character(len=*), intent(in)  datatype,
integer, dimension(6), intent(in), optional  date,
real(dp), intent(in), optional  subsec,
logical, intent(in), optional  append,
logical, intent(in), optional  nohalo 
)

Write 3D data to file.

Parameters
[in]varvalue of the variable
[in]basenamebasename of the file
[in]titletitle of the file
[in]varnamename of the variable
[in]descdescription of the variable
[in]unitunit of the variable
[in]axistypeaxis type (Z/X/Y)
[in]datatypedata type (REAL8/REAL4/default)
[in]dateymdhms of the time
[in]subsecsubsec of the time
[in]appendappend existing (closed) file?
[in]nohaloinclude halo data?

Definition at line 1142 of file scale_fileio.F90.

References gtool_file::filecreate(), fileio_create(), fileio_def_var(), fileio_enddef(), fileio_write_var_3d(), scale_stdio::io_fid_log, scale_stdio::io_l, scale_process::prc_masterrank, scale_process::prc_mpistop(), and gtool_file::rmiss.

1142  use gtool_file, only: &
1143  rmiss
1144  use gtool_file, only: &
1145  filecreate, &
1146  fileaddvariable, &
1147  filesetglobalattribute, &
1148  filewrite
1149  use scale_process, only: &
1150  prc_masterrank, &
1151  prc_mpistop
1152  implicit none
1153 
1154  real(RP), intent(in) :: var(:,:,:)
1155  character(len=*), intent(in) :: basename
1156  character(len=*), intent(in) :: title
1157  character(len=*), intent(in) :: varname
1158  character(len=*), intent(in) :: desc
1159  character(len=*), intent(in) :: unit
1160  character(len=*), intent(in) :: axistype
1161  character(len=*), intent(in) :: datatype
1162 
1163  integer, intent(in), optional :: date(6)
1164  real(DP), intent(in), optional :: subsec
1165  logical, intent(in), optional :: append
1166  logical, intent(in), optional :: nohalo
1167 
1168  integer :: fid, vid
1169  !---------------------------------------------------------------------------
1170 
1171  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Write to file (3D), name : ', trim(varname)
1172 
1173  call fileio_create( fid, & ! [OUT]
1174  basename, title, datatype, date, subsec, append )
1175 
1176  call fileio_def_var( fid, vid, varname, desc, unit, axistype, datatype )
1177 
1178  call fileio_enddef( fid )
1179 
1180  call fileio_write_var_3d( fid, vid, var, varname, axistype, nohalo )
1181 
1182  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
real(dp), parameter, public rmiss
Definition: gtool_file.f90:150
module PROCESS
integer, parameter, public prc_masterrank
master process in each communicator
subroutine, public filecreate(fid, existed, basename, title, source, institution, master, myrank, rankidx, single, time_units, append, comm)
Definition: gtool_file.f90:191
Here is the call graph for this function:

◆ fileio_write_3d_t()

subroutine scale_fileio::fileio_write_3d_t ( real(rp), dimension(:,:,:), intent(in)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  title,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  desc,
character(len=*), intent(in)  unit,
character(len=*), intent(in)  axistype,
character(len=*), intent(in)  datatype,
real(rp), intent(in)  timeintv,
integer, dimension(6), intent(in)  tsince,
logical, intent(in), optional  append,
integer, intent(in), optional  timetarg,
logical, intent(in), optional  nohalo 
)

Write 3D data with time dimension to file.

Parameters
[in]varvalue of the variable
[in]basenamebasename of the file
[in]titletitle of the file
[in]varnamename of the variable
[in]descdescription of the variable
[in]unitunit of the variable
[in]axistypeaxis type (X/Y/Time)
[in]datatypedata type (REAL8/REAL4/default)
[in]timeintvtime interval [sec]
[in]tsincestart time
[in]appendappend existing (closed) file?
[in]timetargtarget timestep (optional)
[in]nohaloinclude halo data?

Definition at line 1201 of file scale_fileio.F90.

References gtool_file_h::file_real4, gtool_file_h::file_real8, gtool_file::filecreate(), fileio_create(), fileio_def_var(), fileio_enddef(), fileio_write_var_3d_t(), scale_stdio::io_fid_log, scale_stdio::io_l, scale_process::prc_masterrank, and scale_process::prc_mpistop().

1201  use gtool_file_h, only: &
1202  file_real8, &
1203  file_real4
1204  use gtool_file, only: &
1205  filecreate, &
1206  fileputaxis, &
1207  fileaddvariable, &
1208  filewrite
1209  use scale_process, only: &
1210  prc_masterrank, &
1211  prc_mpistop
1212  implicit none
1213 
1214  real(RP), intent(in) :: var(:,:,:)
1215  character(len=*), intent(in) :: basename
1216  character(len=*), intent(in) :: title
1217  character(len=*), intent(in) :: varname
1218  character(len=*), intent(in) :: desc
1219  character(len=*), intent(in) :: unit
1220  character(len=*), intent(in) :: axistype
1221  character(len=*), intent(in) :: datatype
1222  real(RP), intent(in) :: timeintv
1223  integer , intent(in) :: tsince(6)
1224 
1225  logical, intent(in), optional :: append
1226  integer, intent(in), optional :: timetarg
1227  logical, intent(in), optional :: nohalo
1228 
1229  integer :: fid, vid
1230  integer :: nsteps
1231 
1232  intrinsic :: size
1233  !---------------------------------------------------------------------------
1234 
1235  if( io_l ) write(io_fid_log,'(1x,3A)') '*** Write to file (3D), name : ', trim(varname), 'with time dimension'
1236 
1237  call fileio_create( fid, & ! [OUT]
1238  basename, title, datatype, tsince, append=append )
1239 
1240  if ( present(timetarg) ) then
1241  nsteps = 1
1242  else
1243  nsteps = size(var,3)
1244  end if
1245  call fileio_def_var( fid, vid, varname, desc, unit, axistype, datatype, timeintv, nsteps )
1246 
1247  call fileio_enddef( fid )
1248 
1249  call fileio_write_var_3d_t( fid, vid, var, varname, axistype, timeintv, timetarg, nohalo )
1250 
1251  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, parameter, public prc_masterrank
master process in each communicator
integer, parameter, public file_real4
subroutine, public filecreate(fid, existed, basename, title, source, institution, master, myrank, rankidx, single, time_units, append, comm)
Definition: gtool_file.f90:191
module FILE I/O HEADER
integer, parameter, public file_real8
Here is the call graph for this function:

◆ fileio_write_4d()

subroutine scale_fileio::fileio_write_4d ( real(rp), dimension(:,:,:,:), intent(in)  var,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  title,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  desc,
character(len=*), intent(in)  unit,
character(len=*), intent(in)  axistype,
character(len=*), intent(in)  datatype,
real(rp), intent(in)  timeintv,
integer, dimension(6), intent(in)  tsince,
logical, intent(in), optional  append,
integer, intent(in), optional  timetarg,
logical, intent(in), optional  nohalo 
)

Write 4D data to file.

Parameters
[in]varvalue of the variable
[in]basenamebasename of the file
[in]titletitle of the file
[in]varnamename of the variable
[in]descdescription of the variable
[in]unitunit of the variable
[in]axistypeaxis type (Z/X/Y/Time)
[in]datatypedata type (REAL8/REAL4/default)
[in]timeintvtime interval [sec]
[in]tsincestart time
[in]appendappend existing (closed) file?
[in]timetargtarget timestep (optional)
[in]nohaloinclude halo data?

Definition at line 1270 of file scale_fileio.F90.

References gtool_file::filecreate(), fileio_create(), fileio_def_var(), fileio_enddef(), fileio_write_var_4d(), scale_stdio::io_fid_log, scale_stdio::io_l, scale_process::prc_mpistop(), and gtool_file::rmiss.

1270  use gtool_file, only: &
1271  rmiss
1272  use gtool_file, only: &
1273  filecreate, &
1274  fileputaxis, &
1275  fileaddvariable, &
1276  filewrite
1277  use scale_process, only: &
1278  prc_mpistop
1279  implicit none
1280 
1281  real(RP), intent(in) :: var(:,:,:,:)
1282  character(len=*), intent(in) :: basename
1283  character(len=*), intent(in) :: title
1284  character(len=*), intent(in) :: varname
1285  character(len=*), intent(in) :: desc
1286  character(len=*), intent(in) :: unit
1287  character(len=*), intent(in) :: axistype
1288  character(len=*), intent(in) :: datatype
1289  real(RP), intent(in) :: timeintv
1290  integer, intent(in) :: tsince(6)
1291 
1292  logical, intent(in), optional :: append
1293  integer, intent(in), optional :: timetarg
1294  logical, intent(in), optional :: nohalo
1295 
1296  integer :: fid, vid
1297  integer :: nsteps
1298 
1299  intrinsic :: size
1300  !---------------------------------------------------------------------------
1301 
1302  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Write to file (4D), name : ', trim(varname)
1303 
1304  call fileio_create( fid, & ! [OUT]
1305  basename, title, datatype, tsince, append=append )
1306 
1307  if ( present(timetarg) ) then
1308  nsteps = 1
1309  else
1310  nsteps = size(var,3)
1311  end if
1312  call fileio_def_var( fid, vid, varname, desc, unit, axistype, datatype, timeintv, nsteps )
1313 
1314  call fileio_enddef( fid )
1315 
1316  call fileio_write_var_4d( fid, vid, var, varname, axistype, timeintv, timetarg, nohalo )
1317 
1318  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
real(dp), parameter, public rmiss
Definition: gtool_file.f90:150
module PROCESS
subroutine, public filecreate(fid, existed, basename, title, source, institution, master, myrank, rankidx, single, time_units, append, comm)
Definition: gtool_file.f90:191
Here is the call graph for this function:

◆ fileio_open()

subroutine, public scale_fileio::fileio_open ( integer, intent(out)  fid,
character(len=*), intent(in)  basename 
)

open a netCDF file for read

Parameters
[out]fidfile ID
[in]basenamebasename of the file

Definition at line 1326 of file scale_fileio.F90.

References gtool_file_h::file_fread, gtool_file::fileopen(), scale_stdio::io_aggregate, scale_process::prc_local_comm_world, scale_process::prc_myrank, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_boundary::atmos_boundary_resume(), mod_atmos_dyn_vars::atmos_dyn_vars_restart_open(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_open(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_open(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_open(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_open(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_open(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_open(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_open(), scale_atmos_refstate::atmos_refstate_read(), mod_atmos_vars::atmos_vars_restart_check(), mod_atmos_vars::atmos_vars_restart_open(), fileio_check_coordinates_name(), fileio_read_1d(), fileio_read_2d(), fileio_read_3d(), fileio_read_4d(), mod_land_vars::land_vars_restart_open(), scale_landuse::landuse_calc_fact(), mod_ocean_vars::ocean_vars_restart_open(), scale_topography::topo_fillhalo(), and mod_urban_vars::urban_vars_restart_open().

1326  use gtool_file_h, only: &
1327  file_fread
1328  use gtool_file, only: &
1329  fileopen
1330  use scale_process, only: &
1331  prc_myrank, &
1333  use mpi, only : mpi_comm_null
1334  implicit none
1335 
1336  integer, intent(out) :: fid
1337  character(len=*), intent(in) :: basename
1338 
1339  integer :: comm
1340  !---------------------------------------------------------------------------
1341 
1342  call prof_rapstart('FILE_O_NetCDF', 2)
1343 
1344  if ( io_aggregate ) then ! user input parameter indicates to do PnetCDF I/O
1345  comm = prc_local_comm_world
1346  else
1347  comm = mpi_comm_null
1348  end if
1349 
1350  call fileopen( fid, & ! [OUT]
1351  basename, & ! [IN]
1352  file_fread, & ! [IN]
1353  comm = comm, & ! [IN]
1354  myrank = prc_myrank ) ! [IN]
1355 
1356  file_closed(fid) = .false.
1357 
1358  call prof_rapend ('FILE_O_NetCDF', 2)
1359 
1360  return
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public prc_local_comm_world
local communicator
module PROCESS
integer, parameter, public file_fread
integer, public prc_myrank
process num in local communicator
module FILE I/O HEADER
subroutine, public fileopen(fid, basename, mode, single, comm, myrank)
Definition: gtool_file.f90:495
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_create()

subroutine, public scale_fileio::fileio_create ( integer, intent(out)  fid,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  title,
character(len=*), intent(in)  datatype,
integer, dimension(6), intent(in), optional  date,
real(dp), intent(in), optional  subsec,
logical, intent(in), optional  append,
logical, intent(in), optional  nozcoord 
)

Create/open a netCDF file.

Parameters
[out]fidfile ID
[in]basenamebasename of the file
[in]titletitle of the file
[in]datatypedata type (REAL8/REAL4/default)
[in]dateymdhms of the time
[in]subsecsubsec of the time
[in]appendswitch whether append existing file or not (default=false)
[in]nozcoordswitch whether include zcoordinate or not (default=false)

Definition at line 1374 of file scale_fileio.F90.

References gtool_file_h::file_real4, gtool_file_h::file_real8, gtool_file::filecreate(), fileio_def_axes(), scale_stdio::h_institute, scale_stdio::h_source, scale_grid_index::ihalo, scale_stdio::io_aggregate, scale_grid_index::jhalo, scale_rm_process::prc_2drank, scale_process::prc_local_comm_world, scale_process::prc_masterrank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_rm_process::prc_periodic_x, scale_rm_process::prc_periodic_y, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_precision::rp, scale_time::time_nowdate, and scale_time::time_nowms.

Referenced by mod_atmos_dyn_vars::atmos_dyn_vars_restart_create(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_create(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_create(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_create(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_create(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_create(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_create(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_create(), mod_atmos_vars::atmos_vars_restart_create(), fileio_write_1d(), fileio_write_2d(), fileio_write_3d(), fileio_write_3d_t(), fileio_write_4d(), mod_land_vars::land_vars_restart_create(), scale_landuse::landuse_write(), mod_ocean_vars::ocean_vars_restart_create(), mod_realinput::parentatomsetup(), and mod_urban_vars::urban_vars_restart_create().

1374  use mpi, only: &
1375  mpi_comm_null
1376  use gtool_file_h, only: &
1377  file_real8, &
1378  file_real4
1379  use gtool_file, only: &
1380  filecreate, &
1381  filesetglobalattribute
1382  use scale_process, only: &
1383  prc_masterrank, &
1384  prc_myrank, &
1385  prc_mpistop, &
1387  use scale_rm_process, only: &
1388  prc_2drank, &
1389  prc_periodic_x, &
1391  use scale_time, only: &
1392  nowdate => time_nowdate, &
1393  nowms => time_nowms
1394  implicit none
1395 
1396  integer, intent(out) :: fid
1397  character(len=*), intent(in) :: basename
1398  character(len=*), intent(in) :: title
1399  character(len=*), intent(in) :: datatype
1400 
1401  integer, intent(in), optional :: date(6)
1402  real(DP), intent(in), optional :: subsec
1403  logical, intent(in), optional :: append
1404  logical, intent(in), optional :: nozcoord
1405 
1406  integer :: rankidx(2)
1407  integer :: dtype
1408  logical :: append_sw
1409  character(len=34) :: tunits
1410  integer :: comm
1411  logical :: fileexisted
1412  character(len=H_SHORT) :: logical_str
1413  !---------------------------------------------------------------------------
1414 
1415  call prof_rapstart('FILE_O_NetCDF', 2)
1416 
1417  rankidx(1) = prc_2drank(prc_myrank,1)
1418  rankidx(2) = prc_2drank(prc_myrank,2)
1419 
1420  ! dtype is used to define the data type of axis variables in file
1421  if ( datatype == 'REAL8' ) then
1422  dtype = file_real8
1423  elseif( datatype == 'REAL4' ) then
1424  dtype = file_real4
1425  else
1426  if ( rp == 8 ) then
1427  dtype = file_real8
1428  elseif( rp == 4 ) then
1429  dtype = file_real4
1430  else
1431  write(*,*) 'xxx unsupported data type. Check!', trim(datatype)
1432  call prc_mpistop
1433  endif
1434  endif
1435 
1436  append_sw = .false.
1437  if ( present(append) ) then
1438  append_sw = append
1439  endif
1440 
1441  ! create a netCDF file if not already existed. Otherwise, open it.
1442  if ( present(date) ) then
1443  call getcftunits( tunits, date )
1444  else
1445  tunits = 'seconds'
1446  end if
1447 
1448  if ( io_aggregate ) then ! user input parameter indicates to do PnetCDF I/O
1449  comm = prc_local_comm_world
1450  else
1451  comm = mpi_comm_null
1452  end if
1453 
1454  call filecreate( fid, & ! [OUT]
1455  fileexisted, & ! [OUT]
1456  basename, & ! [IN]
1457  title, & ! [IN]
1458  h_source, & ! [IN]
1459  h_institute, & ! [IN]
1460  prc_masterrank, & ! [IN]
1461  prc_myrank, & ! [IN]
1462  rankidx, & ! [IN]
1463  time_units = tunits, & ! [IN]
1464  append = append_sw, & ! [IN]
1465  comm = comm ) ! [IN]
1466 
1467  if ( .NOT. fileexisted ) then ! do below only once when file is created
1468  call fileio_def_axes( fid, dtype, xy = nozcoord ) ! [IN]
1469  file_axes_written(fid) = .false. ! indicating axes have not been written yet
1470  if ( present( nozcoord ) ) then
1471  file_nozcoord(fid) = nozcoord
1472  else
1473  file_nozcoord(fid) = .false.
1474  endif
1475 
1476  if ( present( subsec ) ) then
1477  call filesetglobalattribute( fid, "time", (/subsec/) )
1478  else
1479  call filesetglobalattribute( fid, "time", (/nowms/) )
1480  end if
1481  if ( present( date ) ) then
1482  call getcftunits(tunits, date)
1483  else
1484  call getcftunits(tunits, nowdate)
1485  end if
1486  call filesetglobalattribute( fid, "time_units", tunits )
1487  call filesetglobalattribute( fid, "IHALO", (/ihalo/) )
1488  call filesetglobalattribute( fid, "JHALO", (/jhalo/) )
1489  logical_str = "false"
1490  if(prc_periodic_x .AND. .NOT. io_aggregate) logical_str = "true"
1491  call filesetglobalattribute( fid, "PRC_PERIODIC_X", trim(logical_str) )
1492  logical_str = "false"
1493  if(prc_periodic_y .AND. .NOT. io_aggregate) logical_str = "true"
1494  call filesetglobalattribute( fid, "PRC_PERIODIC_Y", trim(logical_str) )
1495 
1496  file_closed(fid) = .false.
1497  endif
1498 
1499  call prof_rapend ('FILE_O_NetCDF', 2)
1500 
1501  return
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public prc_local_comm_world
local communicator
subroutine, public prc_mpistop
Abort MPI.
logical, public prc_periodic_y
periodic condition or not (Y)?
real(dp), public time_nowms
subsecond part of current time [millisec]
Definition: scale_time.F90:66
logical, public prc_periodic_x
periodic condition or not (X)?
integer, public jhalo
of halo cells: y
module TIME
Definition: scale_time.F90:15
module PROCESS
integer, parameter, public prc_masterrank
master process in each communicator
integer, public prc_myrank
process num in local communicator
module RM PROCESS
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, parameter, public file_real4
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
subroutine, public filecreate(fid, existed, basename, title, source, institution, master, myrank, rankidx, single, time_units, append, comm)
Definition: gtool_file.f90:191
module FILE I/O HEADER
integer, parameter, public file_real8
integer, public ihalo
of halo cells: x
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_enddef()

subroutine, public scale_fileio::fileio_enddef ( integer, intent(in)  fid)

Exit netCDF file define mode.

Parameters
[in]fidfile ID

Definition at line 1508 of file scale_fileio.F90.

References gtool_file::fileattachbuffer(), gtool_file::fileenddef(), gtool_file::fileflush(), fileio_write_axes(), scale_stdio::io_aggregate, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by mod_atmos_dyn_vars::atmos_dyn_vars_restart_enddef(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_enddef(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_enddef(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_enddef(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_enddef(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_enddef(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_enddef(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_enddef(), mod_atmos_vars::atmos_vars_restart_enddef(), fileio_write_1d(), fileio_write_2d(), fileio_write_3d(), fileio_write_3d_t(), fileio_write_4d(), mod_land_vars::land_vars_restart_enddef(), scale_landuse::landuse_write(), mod_ocean_vars::ocean_vars_restart_enddef(), mod_realinput::parentatomsetup(), and mod_urban_vars::urban_vars_restart_enddef().

1508  use gtool_file, only: &
1509  fileenddef, &
1510  fileflush, &
1512  implicit none
1513 
1514  integer, intent(in) :: fid
1515  !---------------------------------------------------------------------------
1516 
1517  call prof_rapstart('FILE_O_NetCDF', 2)
1518 
1519  call fileenddef( fid ) ! [IN]
1520 
1521  ! If this enddef is called the first time, write axis variables
1522  if ( .NOT. file_axes_written(fid) ) then
1523  call fileio_write_axes( fid, file_nozcoord(fid) )
1524  if ( io_aggregate ) then
1525  call fileflush( fid )
1526  ! Tell PnetCDF library to use a buffer of size write_buf_amount to
1527  ! aggregate write requests to be post in FILEIO_write_var
1528  call fileattachbuffer( fid, write_buf_amount )
1529  end if
1530  file_axes_written(fid) = .true.
1531  end if
1532 
1533  call prof_rapend ('FILE_O_NetCDF', 2)
1534 
1535  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public fileenddef(fid)
subroutine, public fileattachbuffer(fid, buf_amount)
subroutine, public fileflush(fid)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_flush()

subroutine, public scale_fileio::fileio_flush ( integer, intent(in)  fid)

Flush all pending requests to a netCDF file (PnetCDF only)

Parameters
[in]fidfile ID

Definition at line 1542 of file scale_fileio.F90.

References gtool_file::fileflush(), scale_stdio::io_aggregate, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_boundary::atmos_boundary_update(), mod_atmos_dyn_vars::atmos_dyn_vars_restart_read(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_read(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_read(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_read(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_read(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_read(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_read(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_read(), mod_atmos_vars::atmos_vars_restart_check(), mod_atmos_vars::atmos_vars_restart_read(), fileio_check_coordinates_name(), mod_land_vars::land_vars_restart_read(), scale_landuse::landuse_calc_fact(), mod_ocean_vars::ocean_vars_restart_read(), scale_topography::topo_fillhalo(), and mod_urban_vars::urban_vars_restart_read().

1542  use gtool_file, only: &
1543  fileflush
1544  implicit none
1545 
1546  integer, intent(in) :: fid
1547  !---------------------------------------------------------------------------
1548 
1549  call prof_rapstart('FILE_O_NetCDF', 2)
1550 
1551  if ( io_aggregate ) then
1552  call fileflush( fid ) ! flush all pending read/write requests
1553  end if
1554 
1555  call prof_rapend ('FILE_O_NetCDF', 2)
1556 
1557  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public fileflush(fid)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_close()

subroutine, public scale_fileio::fileio_close ( integer, intent(in)  fid)

Close a netCDF file.

Parameters
[in]fidfile ID

Definition at line 1564 of file scale_fileio.F90.

References gtool_file::fileclose(), gtool_file::filedetachbuffer(), gtool_file::fileflush(), scale_stdio::io_aggregate, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_boundary::atmos_boundary_finalize(), scale_atmos_boundary::atmos_boundary_resume(), mod_atmos_dyn_vars::atmos_dyn_vars_restart_close(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_close(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_close(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_close(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_close(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_close(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_close(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_close(), scale_atmos_refstate::atmos_refstate_read(), mod_atmos_vars::atmos_vars_restart_check(), mod_atmos_vars::atmos_vars_restart_close(), fileio_read_1d(), fileio_read_2d(), fileio_read_3d(), fileio_read_4d(), fileio_write_var_4d(), mod_land_vars::land_vars_restart_close(), scale_landuse::landuse_calc_fact(), scale_landuse::landuse_write(), mod_ocean_vars::ocean_vars_restart_close(), scale_topography::topo_fillhalo(), and mod_urban_vars::urban_vars_restart_close().

1564  use gtool_file, only: &
1565  fileclose, &
1566  fileflush, &
1568  implicit none
1569 
1570  integer, intent(in) :: fid
1571 
1572  !---------------------------------------------------------------------------
1573 
1574  call prof_rapstart('FILE_O_NetCDF', 2)
1575 
1576  if ( .NOT. file_closed(fid) ) then
1577 
1578  if ( io_aggregate ) then
1579  call fileflush( fid ) ! flush all pending read/write requests
1580  if ( write_buf_amount > 0 ) then
1581  call filedetachbuffer( fid ) ! detach PnetCDF aggregation buffer
1582  write_buf_amount = 0 ! reset write request amount
1583  end if
1584  end if
1585 
1586  call fileclose( fid ) ! [IN]
1587 
1588  file_closed(fid) = .true.
1589 
1590  end if
1591 
1592  call prof_rapend ('FILE_O_NetCDF', 2)
1593 
1594  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public fileclose(fid)
subroutine, public filedetachbuffer(fid)
subroutine, public fileflush(fid)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_def_axes()

subroutine scale_fileio::fileio_def_axes ( integer, intent(in)  fid,
integer, intent(in)  dtype,
logical, intent(in), optional  xy 
)

define axis variables in the file

Definition at line 1603 of file scale_fileio.F90.

References gtool_file::filedefassociatedcoordinates(), gtool_file::filedefaxis(), gtool_file::filesettattr(), scale_grid_index::ia, scale_grid_index::iag, scale_grid_index::imaxb, scale_stdio::io_aggregate, scale_grid_index::ja, scale_grid_index::jag, scale_grid_index::jmaxb, scale_grid_index::ka, scale_grid_index::kmax, scale_land_grid_index::lkmax, scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, and scale_urban_grid_index::ukmax.

Referenced by fileio_create().

1603  use gtool_file, only: &
1604  filedefaxis, &
1605  filesettattr, &
1607  use scale_rm_process, only: &
1608  prc_num_x, &
1609  prc_num_y
1610  implicit none
1611 
1612  integer, intent(in) :: fid
1613  integer, intent(in) :: dtype
1614  logical, intent(in), optional :: xy
1615 
1616  character(len=2) :: AXIS_name(3)
1617  logical :: xy_
1618  !---------------------------------------------------------------------------
1619 
1620  if ( present(xy) ) then
1621  xy_ = xy
1622  else
1623  xy_ = .false.
1624  end if
1625 
1626  if ( .NOT. xy_ ) then
1627  call filedefaxis( fid, 'z', 'Z', 'm', 'z', dtype, kmax )
1628  end if
1629  if ( .NOT. io_aggregate ) then
1630  call filedefaxis( fid, 'x', 'X', 'm', 'x', dtype, imaxb )
1631  call filedefaxis( fid, 'y', 'Y', 'm', 'y', dtype, jmaxb )
1632  else
1633  call filedefaxis( fid, 'x', 'X', 'm', 'x', dtype, iag )
1634  call filedefaxis( fid, 'y', 'Y', 'm', 'y', dtype, jag )
1635  end if
1636 
1637  if ( .NOT. xy_ ) then
1638  call filedefaxis( fid, 'zh', 'Z (half level)', 'm', 'zh', dtype, kmax )
1639  end if
1640  if ( .NOT. io_aggregate ) then
1641  call filedefaxis( fid, 'xh', 'X (half level)', 'm', 'xh', dtype, imaxb )
1642  call filedefaxis( fid, 'yh', 'Y (half level)', 'm', 'yh', dtype, jmaxb )
1643  else
1644  call filedefaxis( fid, 'xh', 'X (half level)', 'm', 'xh', dtype, iag )
1645  call filedefaxis( fid, 'yh', 'Y (half level)', 'm', 'yh', dtype, jag )
1646  end if
1647 
1648  if ( .NOT. xy_ ) then
1649  call filedefaxis( fid, 'lz', 'LZ', 'm', 'lz', dtype, lkmax )
1650  call filedefaxis( fid, 'lzh', 'LZ (half level)', 'm', 'lzh', dtype, lkmax )
1651  call filedefaxis( fid, 'uz', 'UZ', 'm', 'uz', dtype, ukmax )
1652  call filedefaxis( fid, 'uzh', 'UZ (half level)', 'm', 'uzh', dtype, ukmax )
1653  end if
1654 
1655  if ( .NOT. xy_ ) then
1656  call filedefaxis( fid, 'CZ', 'Atmos Grid Center Position Z', 'm', 'CZ', dtype, ka )
1657  end if
1658  if ( .NOT. io_aggregate ) then
1659  call filedefaxis( fid, 'CX', 'Atmos Grid Center Position X', 'm', 'CX', dtype, ia )
1660  call filedefaxis( fid, 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', dtype, ja )
1661  else
1662  call filedefaxis( fid, 'CX', 'Atmos Grid Center Position X', 'm', 'CX', dtype, iag )
1663  call filedefaxis( fid, 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', dtype, jag )
1664  end if
1665 
1666  if ( .NOT. xy_ ) then
1667  call filedefaxis( fid, 'FZ', 'Atmos Grid Face Position Z', 'm', 'FZ', dtype, ka+1 )
1668  end if
1669  if ( .NOT. io_aggregate ) then
1670  call filedefaxis( fid, 'FX', 'Atmos Grid Face Position X', 'm', 'FX', dtype, ia+1 )
1671  call filedefaxis( fid, 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', dtype, ja+1 )
1672  else
1673  call filedefaxis( fid, 'FX', 'Atmos Grid Face Position X', 'm', 'FX', dtype, iag+1 )
1674  call filedefaxis( fid, 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', dtype, jag+1 )
1675  end if
1676 
1677  if ( .NOT. xy_ ) then
1678  call filedefaxis( fid, 'CDZ', 'Grid Cell length Z', 'm', 'CZ', dtype, ka )
1679  end if
1680  if ( .NOT. io_aggregate ) then
1681  call filedefaxis( fid, 'CDX', 'Grid Cell length X', 'm', 'CX', dtype, ia )
1682  call filedefaxis( fid, 'CDY', 'Grid Cell length Y', 'm', 'CY', dtype, ja )
1683  else
1684  call filedefaxis( fid, 'CDX', 'Grid Cell length X', 'm', 'CX', dtype, iag )
1685  call filedefaxis( fid, 'CDY', 'Grid Cell length Y', 'm', 'CY', dtype, jag )
1686  end if
1687 
1688  if ( .NOT. xy_ ) then
1689  call filedefaxis( fid, 'FDZ', 'Grid distance Z', 'm', 'FDZ', dtype, ka-1 )
1690  end if
1691  if ( .NOT. io_aggregate ) then
1692  call filedefaxis( fid, 'FDX', 'Grid distance X', 'm', 'FDX', dtype, ia-1 )
1693  call filedefaxis( fid, 'FDY', 'Grid distance Y', 'm', 'FDY', dtype, ja-1 )
1694  else
1695  call filedefaxis( fid, 'FDX', 'Grid distance X', 'm', 'FDX', dtype, iag-1 )
1696  call filedefaxis( fid, 'FDY', 'Grid distance Y', 'm', 'FDY', dtype, jag-1 )
1697  end if
1698 
1699  if ( .NOT. xy_ ) then
1700  call filedefaxis( fid, 'LCZ', 'Land Grid Center Position Z', 'm', 'LCZ', dtype, lkmax )
1701  call filedefaxis( fid, 'LFZ', 'Land Grid Face Position Z', 'm', 'LFZ', dtype, lkmax+1 )
1702  call filedefaxis( fid, 'LCDZ', 'Land Grid Cell length Z', 'm', 'LCZ', dtype, lkmax )
1703 
1704  call filedefaxis( fid, 'UCZ', 'Urban Grid Center Position Z', 'm', 'UCZ', dtype, ukmax )
1705  call filedefaxis( fid, 'UFZ', 'Urban Grid Face Position Z', 'm', 'UFZ', dtype, ukmax+1 )
1706  call filedefaxis( fid, 'UCDZ', 'Urban Grid Cell length Z', 'm', 'UCZ', dtype, ukmax )
1707  end if
1708 
1709  if ( .NOT. xy_ ) then
1710  call filedefaxis( fid, 'CBFZ', 'Boundary factor Center Z', '1', 'CZ', dtype, ka )
1711  end if
1712  if ( .NOT. io_aggregate ) then
1713  call filedefaxis( fid, 'CBFX', 'Boundary factor Center X', '1', 'CX', dtype, ia )
1714  call filedefaxis( fid, 'CBFY', 'Boundary factor Center Y', '1', 'CY', dtype, ja )
1715  else
1716  call filedefaxis( fid, 'CBFX', 'Boundary factor Center X', '1', 'CX', dtype, iag )
1717  call filedefaxis( fid, 'CBFY', 'Boundary factor Center Y', '1', 'CY', dtype, jag )
1718  end if
1719 
1720  if ( .NOT. xy_ ) then
1721  call filedefaxis( fid, 'FBFZ', 'Boundary factor Face Z', '1', 'CZ', dtype, ka )
1722  end if
1723  if ( .NOT. io_aggregate ) then
1724  call filedefaxis( fid, 'FBFX', 'Boundary factor Face X', '1', 'CX', dtype, ia )
1725  call filedefaxis( fid, 'FBFY', 'Boundary factor Face Y', '1', 'CY', dtype, ja )
1726  else
1727  call filedefaxis( fid, 'FBFX', 'Boundary factor Face X', '1', 'CX', dtype, iag )
1728  call filedefaxis( fid, 'FBFY', 'Boundary factor Face Y', '1', 'CY', dtype, jag )
1729  end if
1730 
1731  ! TODO: skip 8 axes below when IO_AGGREGATE is true, as all axes are now global
1732  call filedefaxis( fid, 'CXG', 'Grid Center Position X (global)', 'm', 'CXG', dtype, iag )
1733  call filedefaxis( fid, 'CYG', 'Grid Center Position Y (global)', 'm', 'CYG', dtype, jag )
1734  call filedefaxis( fid, 'FXG', 'Grid Face Position X (global)', 'm', 'FXG', dtype, iag+1 )
1735  call filedefaxis( fid, 'FYG', 'Grid Face Position Y (global)', 'm', 'FYG', dtype, jag+1 )
1736 
1737  call filedefaxis( fid, 'CBFXG', 'Boundary factor Center X (global)', '1', 'CXG', dtype, iag )
1738  call filedefaxis( fid, 'CBFYG', 'Boundary factor Center Y (global)', '1', 'CYG', dtype, jag )
1739  call filedefaxis( fid, 'FBFXG', 'Boundary factor Face X (global)', '1', 'CXG', dtype, iag )
1740  call filedefaxis( fid, 'FBFYG', 'Boundary factor Face Y (global)', '1', 'CYG', dtype, jag )
1741 
1742  ! associate coordinates
1743  axis_name(1:2) = (/'x ','y '/)
1744  call filedefassociatedcoordinates( fid, 'lon' , 'longitude', &
1745  'degrees_east' , axis_name(1:2), dtype )
1746  axis_name(1:2) = (/'xh','y '/)
1747  call filedefassociatedcoordinates( fid, 'lon_uy', 'longitude (half level uy)', &
1748  'degrees_east' , axis_name(1:2), dtype )
1749  axis_name(1:2) = (/'x ','yh'/)
1750  call filedefassociatedcoordinates( fid, 'lon_xv', 'longitude (half level xv)', &
1751  'degrees_east' , axis_name(1:2), dtype )
1752  axis_name(1:2) = (/'xh','yh'/)
1753  call filedefassociatedcoordinates( fid, 'lon_uv', 'longitude (half level uv)', &
1754  'degrees_east' , axis_name(1:2), dtype )
1755  axis_name(1:2) = (/'x ','y '/)
1756  call filedefassociatedcoordinates( fid, 'lat' , 'latitude', &
1757  'degrees_north', axis_name(1:2), dtype )
1758  axis_name(1:2) = (/'xh','y '/)
1759  call filedefassociatedcoordinates( fid, 'lat_uy', 'latitude (half level uy)', &
1760  'degrees_north', axis_name(1:2), dtype )
1761  axis_name(1:2) = (/'x ','yh'/)
1762  call filedefassociatedcoordinates( fid, 'lat_xv', 'latitude (half level xv)', &
1763  'degrees_north', axis_name(1:2), dtype )
1764  axis_name(1:2) = (/'xh','yh'/)
1765  call filedefassociatedcoordinates( fid, 'lat_uv', 'latitude (half level uv)', &
1766  'degrees_north', axis_name(1:2), dtype )
1767 
1768  if ( .NOT. xy_ ) then
1769  axis_name = (/'z', 'x', 'y'/)
1770  call filedefassociatedcoordinates( fid, 'height', 'height above ground level', &
1771  'm', axis_name(1:3), dtype )
1772  axis_name = (/'zh', 'x ', 'y '/)
1773  call filedefassociatedcoordinates( fid, 'height_wxy', 'height above ground level (half level wxy)', &
1774  'm', axis_name(1:3), dtype )
1775  end if
1776 
1777  ! attributes
1778  if ( .NOT. xy_ ) then
1779  call filesettattr( fid, 'lz', 'positive', 'down' )
1780  call filesettattr( fid, 'lzh', 'positive', 'down' )
1781  call filesettattr( fid, 'uz', 'positive', 'down' )
1782  call filesettattr( fid, 'uzh', 'positive', 'down' )
1783  call filesettattr( fid, 'LCZ', 'positive', 'down' )
1784  call filesettattr( fid, 'LFZ', 'positive', 'down' )
1785  call filesettattr( fid, 'UCZ', 'positive', 'down' )
1786  call filesettattr( fid, 'UFZ', 'positive', 'down' )
1787  end if
1788 
1789  return
integer, public prc_num_x
x length of 2D processor topology
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public imaxb
integer, public jmaxb
integer, public prc_num_y
y length of 2D processor topology
integer, public jag
of computational grids
subroutine, public filedefassociatedcoordinates(fid, name, desc, units, dim_names, dtype)
Definition: gtool_file.f90:903
subroutine, public filesettattr(fid, vname, key, val)
integer, public iag
of computational grids
module RM PROCESS
subroutine, public filedefaxis(fid, name, desc, units, dim_name, dtype, dim_size)
Definition: gtool_file.f90:584
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_write_axes()

subroutine scale_fileio::fileio_write_axes ( integer, intent(in)  fid,
logical, intent(in), optional  xy 
)

write axis to the file

Definition at line 1797 of file scale_fileio.F90.

References scale_grid::grid_cbfx, scale_grid::grid_cbfxg, scale_grid::grid_cbfy, scale_grid::grid_cbfyg, scale_grid::grid_cbfz, scale_grid::grid_cdx, scale_grid::grid_cdxg, scale_grid::grid_cdy, scale_grid::grid_cdyg, scale_grid::grid_cdz, scale_grid::grid_cx, scale_grid::grid_cxg, scale_grid::grid_cy, scale_grid::grid_cyg, scale_grid::grid_cz, scale_grid::grid_fbfx, scale_grid::grid_fbfxg, scale_grid::grid_fbfy, scale_grid::grid_fbfyg, scale_grid::grid_fbfz, scale_grid::grid_fdx, scale_grid::grid_fdxg, scale_grid::grid_fdy, scale_grid::grid_fdyg, scale_grid::grid_fdz, scale_grid::grid_fx, scale_grid::grid_fxg, scale_grid::grid_fy, scale_grid::grid_fyg, scale_grid::grid_fz, scale_land_grid::grid_lcdz, scale_land_grid::grid_lcz, scale_land_grid::grid_lfz, scale_urban_grid::grid_ucdz, scale_urban_grid::grid_ucz, scale_urban_grid::grid_ufz, scale_stdio::io_aggregate, scale_grid_index::isga, scale_grid_index::jsga, scale_grid_index::ke, scale_grid_index::ks, scale_land_grid_index::lke, scale_land_grid_index::lks, scale_rm_process::prc_2drank, scale_process::prc_myrank, scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_urban_grid_index::uke, and scale_urban_grid_index::uks.

Referenced by fileio_enddef().

1797  use gtool_file, only: &
1798  filewriteaxis, &
1799  filewriteassociatedcoordinates
1800  use scale_process, only: &
1801  prc_myrank
1802  use scale_rm_process, only: &
1803  prc_2drank, &
1804  prc_num_x, &
1805  prc_num_y
1806  use scale_grid, only: &
1807  grid_cz, &
1808  grid_cx, &
1809  grid_cy, &
1810  grid_fz, &
1811  grid_fx, &
1812  grid_fy, &
1813  grid_cdz, &
1814  grid_cdx, &
1815  grid_cdy, &
1816  grid_fdz, &
1817  grid_fdx, &
1818  grid_fdy, &
1819  grid_cbfz, &
1820  grid_cbfx, &
1821  grid_cbfy, &
1822  grid_fbfz, &
1823  grid_fbfx, &
1824  grid_fbfy, &
1825  grid_cxg, &
1826  grid_cyg, &
1827  grid_fxg, &
1828  grid_fyg, &
1829  grid_cdxg, &
1830  grid_cdyg, &
1831  grid_fdxg, &
1832  grid_fdyg, &
1833  grid_cbfxg, &
1834  grid_cbfyg, &
1835  grid_fbfxg, &
1836  grid_fbfyg
1837  use scale_land_grid, only: &
1838  grid_lcz, &
1839  grid_lfz, &
1840  grid_lcdz
1841  use scale_urban_grid, only: &
1842  grid_ucz, &
1843  grid_ufz, &
1844  grid_ucdz
1845  implicit none
1846 
1847  integer, intent(in) :: fid
1848  logical, intent(in), optional :: xy
1849 
1850  logical :: xy_
1851  logical :: put_x
1852  logical :: put_y
1853  logical :: put_z
1854 
1855  integer :: rankidx(2)
1856  integer :: start(3)
1857  !---------------------------------------------------------------------------
1858 
1859  if ( present(xy) ) then
1860  xy_ = xy
1861  else
1862  xy_ = .false.
1863  end if
1864 
1865  if ( io_aggregate ) then
1866  rankidx(1) = prc_2drank(prc_myrank,1)
1867  rankidx(2) = prc_2drank(prc_myrank,2)
1868 
1869  ! For parallel I/O, not all variables are written by all processes.
1870  ! 1. Let PRC_myrank 0 writes all z axes
1871  ! 2. Let processes (rankidx(2) == 0) write x axes (south-most processes)
1872  ! rankidx(1) == 0 writes west HALO
1873  ! rankidx(1) == PRC_NUM_X-1 writes east HALO
1874  ! others writes without HALO
1875  ! 3. Let processes (rankidx(1) == 0) write y axes (west-most processes)
1876  ! rankidx(1) == 0 writes south HALO
1877  ! rankidx(1) == PRC_NUM_Y-1 writes north HALO
1878  ! others writes without HALO
1879 
1880  put_z = ( .NOT. xy_ ) .AND. ( prc_myrank == 0 ) ! only master process output the vertical coordinates
1881  put_x = ( rankidx(2) == 0 ) ! only south most row processes write x coordinates
1882  put_y = ( rankidx(1) == 0 ) ! only west most column processes write y coordinates
1883  else
1884  put_z = ( .NOT. xy_ )
1885  put_x = .true.
1886  put_y = .true.
1887 
1888  start(:) = 1
1889  end if
1890 
1891  if ( put_z ) then
1892  if( io_aggregate ) start(1) = 1
1893  call filewriteaxis( fid, 'z', grid_cz(ks:ke), start )
1894  call filewriteaxis( fid, 'zh', grid_fz(ks:ke), start )
1895  call filewriteaxis( fid, 'lz', grid_lcz(lks:lke), start )
1896  call filewriteaxis( fid, 'lzh', grid_lfz(lks:lke), start )
1897  call filewriteaxis( fid, 'uz', grid_ucz(uks:uke), start )
1898  call filewriteaxis( fid, 'uzh', grid_ufz(uks:uke), start )
1899  end if
1900  if ( put_x ) then
1901  if( io_aggregate ) start(1) = isga
1902  call filewriteaxis( fid, 'x', grid_cx(xsb:xeb), start )
1903  call filewriteaxis( fid, 'xh', grid_fx(xsb:xeb), start )
1904  end if
1905  if ( put_y ) then
1906  if( io_aggregate ) start(1) = jsga
1907  call filewriteaxis( fid, 'y', grid_cy(ysb:yeb), start )
1908  call filewriteaxis( fid, 'yh', grid_fy(ysb:yeb), start )
1909  end if
1910 
1911  if ( put_z ) then
1912  if( io_aggregate ) start(1) = 1
1913  call filewriteaxis( fid, 'CZ', grid_cz, start )
1914  call filewriteaxis( fid, 'FZ', grid_fz, start )
1915  call filewriteaxis( fid, 'CDZ', grid_cdz, start )
1916  call filewriteaxis( fid, 'FDZ', grid_fdz, start )
1917 
1918  call filewriteaxis( fid, 'LCZ', grid_lcz, start )
1919  call filewriteaxis( fid, 'LFZ', grid_lfz, start )
1920  call filewriteaxis( fid, 'LCDZ', grid_lcdz, start )
1921 
1922  call filewriteaxis( fid, 'UCZ', grid_ucz, start )
1923  call filewriteaxis( fid, 'UFZ', grid_ufz, start )
1924  call filewriteaxis( fid, 'UCDZ', grid_ucdz, start )
1925 
1926  call filewriteaxis( fid, 'CBFZ', grid_cbfz, start )
1927  call filewriteaxis( fid, 'FBFZ', grid_fbfz, start )
1928  end if
1929 
1930  if ( io_aggregate ) then
1931  if ( prc_myrank == 0 ) then
1932  start(1) = 1
1933  call filewriteaxis( fid, 'CX', grid_cxg, start )
1934  call filewriteaxis( fid, 'CY', grid_cyg, start )
1935  call filewriteaxis( fid, 'FX', grid_fxg, start )
1936  call filewriteaxis( fid, 'FY', grid_fyg, start )
1937 
1938  call filewriteaxis( fid, 'CDX', grid_cdxg, start )
1939  call filewriteaxis( fid, 'CDY', grid_cdyg, start )
1940  call filewriteaxis( fid, 'FDX', grid_fdxg, start )
1941  call filewriteaxis( fid, 'FDY', grid_fdyg, start )
1942 
1943  call filewriteaxis( fid, 'CBFX', grid_cbfxg, start )
1944  call filewriteaxis( fid, 'CBFY', grid_cbfyg, start )
1945  call filewriteaxis( fid, 'FBFX', grid_fbfxg, start )
1946  call filewriteaxis( fid, 'FBFY', grid_fbfyg, start )
1947 
1948  call filewriteaxis( fid, 'CXG', grid_cxg, start )
1949  call filewriteaxis( fid, 'CYG', grid_cyg, start )
1950  call filewriteaxis( fid, 'FXG', grid_fxg, start )
1951  call filewriteaxis( fid, 'FYG', grid_fyg, start )
1952 
1953  call filewriteaxis( fid, 'CBFXG', grid_cbfxg, start )
1954  call filewriteaxis( fid, 'CBFYG', grid_cbfyg, start )
1955  call filewriteaxis( fid, 'FBFXG', grid_fbfxg, start )
1956  call filewriteaxis( fid, 'FBFYG', grid_fbfyg, start )
1957  end if
1958  else
1959  call filewriteaxis( fid, 'CX', grid_cx )
1960  call filewriteaxis( fid, 'CY', grid_cy )
1961  call filewriteaxis( fid, 'FX', grid_fx )
1962  call filewriteaxis( fid, 'FY', grid_fy )
1963 
1964  call filewriteaxis( fid, 'CDX', grid_cdx )
1965  call filewriteaxis( fid, 'CDY', grid_cdy )
1966  call filewriteaxis( fid, 'FDX', grid_fdx )
1967  call filewriteaxis( fid, 'FDY', grid_fdy )
1968 
1969  call filewriteaxis( fid, 'CBFX', grid_cbfx )
1970  call filewriteaxis( fid, 'CBFY', grid_cbfy )
1971  call filewriteaxis( fid, 'FBFX', grid_fbfx )
1972  call filewriteaxis( fid, 'FBFY', grid_fbfy )
1973 
1974  call filewriteaxis( fid, 'CXG', grid_cxg )
1975  call filewriteaxis( fid, 'CYG', grid_cyg )
1976  call filewriteaxis( fid, 'FXG', grid_fxg )
1977  call filewriteaxis( fid, 'FYG', grid_fyg )
1978 
1979  call filewriteaxis( fid, 'CBFXG', grid_cbfxg )
1980  call filewriteaxis( fid, 'CBFYG', grid_cbfyg )
1981  call filewriteaxis( fid, 'FBFXG', grid_fbfxg )
1982  call filewriteaxis( fid, 'FBFYG', grid_fbfyg )
1983  end if
1984 
1985  ! associate coordinates
1986 
1987  if ( io_aggregate ) then
1988  start(1) = isga
1989  start(2) = jsga
1990  end if
1991  call filewriteassociatedcoordinates( fid, 'lon' , axis_lon(xsb:xeb,ysb:yeb), start )
1992  call filewriteassociatedcoordinates( fid, 'lon_uy', axis_lonx(xsb:xeb,ysb:yeb), start )
1993  call filewriteassociatedcoordinates( fid, 'lon_xv', axis_lony(xsb:xeb,ysb:yeb), start )
1994  call filewriteassociatedcoordinates( fid, 'lon_uv', axis_lonxy(xsb:xeb,ysb:yeb), start )
1995  call filewriteassociatedcoordinates( fid, 'lat' , axis_lat(xsb:xeb,ysb:yeb), start )
1996  call filewriteassociatedcoordinates( fid, 'lat_uy', axis_latx(xsb:xeb,ysb:yeb), start )
1997  call filewriteassociatedcoordinates( fid, 'lat_xv', axis_laty(xsb:xeb,ysb:yeb), start )
1998  call filewriteassociatedcoordinates( fid, 'lat_uv', axis_latxy(xsb:xeb,ysb:yeb), start )
1999 
2000  if ( .NOT. xy_ ) then
2001  if ( io_aggregate ) then
2002  start(1) = 1
2003  start(2) = isga
2004  start(3) = jsga
2005  end if
2006  call filewriteassociatedcoordinates( fid, 'height', axis_hzxy(ks:ke,xsb:xeb,ysb:yeb), start )
2007  call filewriteassociatedcoordinates( fid, 'height_wxy', axis_hwxy(ks:ke,xsb:xeb,ysb:yeb), start )
2008  end if
2009 
2010  return
integer, public prc_num_x
x length of 2D processor topology
module GTOOL_FILE
Definition: gtool_file.f90:17
real(rp), dimension(:), allocatable, public grid_cyg
center coordinate [m]: y, global
real(rp), dimension(:), allocatable, public grid_cbfyg
center buffer factor (0-1): y, global
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]
real(rp), dimension(:), allocatable, public grid_cdyg
center coordinate [m]: y, global
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fxg
face coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fbfy
face buffer factor (0-1): y
module GRID (cartesian) for land
real(rp), dimension(:), allocatable, public grid_cdxg
center coordinate [m]: x, global
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 prc_num_y
y length of 2D processor topology
real(rp), dimension(:), allocatable, public grid_ucdz
z-length of control volume [m]
integer, public jsga
start point of the full domain: cy, global
module GRID (cartesian) for urban
real(rp), dimension(:), allocatable, public grid_fbfx
face buffer factor (0-1): x
real(rp), dimension(:), allocatable, public grid_fdz
z-length of grid(k+1) to grid(k) [m]
real(rp), dimension(:), allocatable, public grid_ufz
face coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fbfxg
face buffer factor (0-1): x, global
real(rp), dimension(:), allocatable, public grid_fdyg
center coordinate [m]: y, global
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
real(rp), dimension(:), allocatable, public grid_cbfx
center buffer factor (0-1): x
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
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)
module RM PROCESS
integer, public isga
start point of the full domain: cx, global
real(rp), dimension(:), allocatable, public grid_lcdz
z-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_fdxg
center coordinate [m]: x, global
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
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
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
real(rp), dimension(:), allocatable, public grid_cbfxg
center buffer factor (0-1): x, global
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local
Here is the caller graph for this function:

◆ fileio_def_var()

subroutine, public scale_fileio::fileio_def_var ( integer, intent(in)  fid,
integer, intent(out)  vid,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  desc,
character(len=*), intent(in)  unit,
character(len=*), intent(in)  axistype,
character(len=*), intent(in)  datatype,
real(rp), intent(in), optional  timeintv,
integer, intent(in), optional  nsteps 
)

Define a variable to file.

Parameters
[in]fidfile ID
[out]vidvariable ID
[in]varnamename of the variable
[in]descdescription of the variable
[in]unitunit of the variable
[in]axistypeaxis type (Z/X/Y)
[in]datatypedata type (REAL8/REAL4/default)
[in]timeintvtime interval [sec]
[in]nstepsnumber of time steps

Definition at line 2025 of file scale_fileio.F90.

References gtool_file_h::file_real4, gtool_file_h::file_real8, gtool_file::filedefinevariable(), scale_grid_index::ia, scale_grid_index::ja, scale_grid_index::ka, scale_land_grid_index::lkmax, scale_process::prc_mpistop(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_precision::rp, and scale_urban_grid_index::ukmax.

Referenced by mod_atmos_dyn_vars::atmos_dyn_vars_restart_def_var(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_def_var(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_def_var(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_def_var(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_def_var(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_def_var(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_def_var(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_def_var(), mod_atmos_vars::atmos_vars_restart_def_var(), fileio_write_1d(), fileio_write_2d(), fileio_write_3d(), fileio_write_3d_t(), fileio_write_4d(), mod_land_vars::land_vars_restart_def_var(), scale_landuse::landuse_write(), mod_ocean_vars::ocean_vars_restart_def_var(), mod_realinput::parentatomsetup(), and mod_urban_vars::urban_vars_restart_def_var().

2025  use gtool_file_h, only: &
2026  file_real8, &
2027  file_real4
2028  use gtool_file, only: &
2030  use scale_process, only: &
2031  prc_mpistop
2032  implicit none
2033 
2034  integer, intent(in) :: fid
2035  integer, intent(out) :: vid
2036  character(len=*), intent(in) :: varname
2037  character(len=*), intent(in) :: desc
2038  character(len=*), intent(in) :: unit
2039  character(len=*), intent(in) :: axistype
2040  character(len=*), intent(in) :: datatype
2041 
2042  real(RP), intent(in), optional :: timeintv
2043  integer, intent(in), optional :: nsteps
2044 
2045  integer :: dtype, ndims, elm_size
2046  character(len=2) :: dims(3)
2047  real(DP) :: time_interval
2048  !---------------------------------------------------------------------------
2049 
2050  call prof_rapstart('FILE_O_NetCDF', 2)
2051 
2052  if ( datatype == 'REAL8' ) then
2053  dtype = file_real8
2054  elm_size = 4
2055  elseif( datatype == 'REAL4' ) then
2056  dtype = file_real4
2057  elm_size = 8
2058  else
2059  if ( rp == 8 ) then
2060  dtype = file_real8
2061  elm_size = 8
2062  elseif( rp == 4 ) then
2063  dtype = file_real4
2064  elm_size = 4
2065  else
2066  write(*,*) 'xxx unsupported data type. Check!', trim(datatype), ' item:',trim(varname)
2067  call prc_mpistop
2068  endif
2069  endif
2070 
2071  if ( axistype == 'Z' ) then ! 1D variable
2072  ndims = 1
2073  dims(1) = 'z'
2074  write_buf_amount = write_buf_amount + ka * elm_size
2075  elseif( axistype == 'X' ) then
2076  ndims = 1
2077  dims(1) = 'x'
2078  write_buf_amount = write_buf_amount + ia * elm_size
2079  elseif( axistype == 'Y' ) then
2080  ndims = 1
2081  dims(1) = 'y'
2082  write_buf_amount = write_buf_amount + ja * elm_size
2083  elseif( axistype == 'XY' ) then ! 2D variable
2084  ndims = 2
2085  dims(1) = 'x'
2086  dims(2) = 'y'
2087  write_buf_amount = write_buf_amount + ia * ja * elm_size
2088  elseif( axistype == 'UY' ) then
2089  ndims = 2
2090  dims(1) = 'xh'
2091  dims(2) = 'y'
2092  write_buf_amount = write_buf_amount + ia * ja * elm_size
2093  elseif( axistype == 'XV' ) then
2094  ndims = 2
2095  dims(1) = 'x'
2096  dims(2) = 'yh'
2097  write_buf_amount = write_buf_amount + ia * ja * elm_size
2098  elseif( axistype == 'UV' ) then
2099  ndims = 2
2100  dims(1) = 'xh'
2101  dims(2) = 'yh'
2102  write_buf_amount = write_buf_amount + ia * ja * elm_size
2103  elseif( axistype == 'ZX' ) then
2104  ndims = 2
2105  dims(1) = 'z'
2106  dims(2) = 'x'
2107  write_buf_amount = write_buf_amount + ka * ia * elm_size
2108  elseif( axistype == 'ZXY' ) then ! 3D variable
2109  ndims = 3
2110  dims = (/'z','x','y'/)
2111  write_buf_amount = write_buf_amount + ka * ia * ja * elm_size
2112  elseif( axistype == 'ZHXY' ) then
2113  ndims = 3
2114  dims = (/'zh','x ','y '/)
2115  write_buf_amount = write_buf_amount + ka * ia * ja * elm_size
2116  elseif( axistype == 'ZXHY' ) then
2117  ndims = 3
2118  dims = (/'z ','xh','y '/)
2119  write_buf_amount = write_buf_amount + ka * ia * ja * elm_size
2120  elseif( axistype == 'ZXYH' ) then
2121  ndims = 3
2122  dims = (/'z ','x ','yh'/)
2123  write_buf_amount = write_buf_amount + ka * ia * ja * elm_size
2124  elseif( axistype == 'Land' ) then
2125  ndims = 3
2126  dims = (/'lz','x ','y '/)
2127  write_buf_amount = write_buf_amount + lkmax * ia * ja * elm_size
2128  elseif( axistype == 'Urban' ) then
2129  ndims = 3
2130  dims = (/'uz','x ','y '/)
2131  write_buf_amount = write_buf_amount + ukmax * ia * ja * elm_size
2132  elseif( axistype == 'XYT' ) then ! 3D variable with time dimension
2133  ndims = 2
2134  dims(1) = 'x'
2135  dims(2) = 'y'
2136  if ( present(nsteps) ) then
2137  write_buf_amount = write_buf_amount + ia * ja * elm_size * nsteps
2138  else
2139  write_buf_amount = write_buf_amount + ia * ja * elm_size
2140  end if
2141  elseif( axistype == 'ZXYT' ) then ! 4D variable
2142  ndims = 3
2143  dims = (/'z','x','y'/)
2144  if ( present(nsteps) ) then
2145  write_buf_amount = write_buf_amount + ka * ia * ja * elm_size * nsteps
2146  else
2147  write_buf_amount = write_buf_amount + ka * ia * ja * elm_size
2148  end if
2149  else
2150  write(*,*) 'xxx [FILEIO_def_var] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
2151  call prc_mpistop
2152  endif
2153 
2154  if ( present(timeintv) ) then ! 3D/4D variable with time dimension
2155  time_interval = timeintv
2156  call filedefinevariable( fid, vid, varname, desc, unit, ndims, dims, dtype, & ! [IN]
2157  tint=time_interval ) ! [IN]
2158  else
2159  call filedefinevariable( fid, vid, varname, desc, unit, ndims, dims, dtype ) ! [IN]
2160  endif
2161 
2162  call prof_rapend ('FILE_O_NetCDF', 2)
2163 
2164  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
subroutine, public filedefinevariable(fid, vid, varname, desc, units, ndims, dims, dtype, tint, tavg)
module PROCESS
integer, parameter, public file_real4
module FILE I/O HEADER
integer, parameter, public file_real8
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_write_var_1d()

subroutine scale_fileio::fileio_write_var_1d ( integer, intent(in)  fid,
integer, intent(in)  vid,
real(rp), dimension(:), intent(in)  var,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype 
)

Write 1D data to file.

Parameters
[in]fidfile ID
[in]vidnetCDF variable ID
[in]varvalue of the variable
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y)

Definition at line 2175 of file scale_fileio.F90.

References scale_grid_index::ieb, scale_stdio::io_aggregate, scale_grid_index::isb, scale_grid_index::isga, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::jsga, scale_grid_index::ke, scale_grid_index::ks, scale_rm_process::prc_2drank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_time::time_nowdaysec.

Referenced by fileio_write_1d().

2175  use gtool_file, only: &
2176  filewrite
2177  use scale_process, only: &
2178  prc_myrank, &
2179  prc_mpistop
2180  use scale_rm_process, only: &
2181  prc_2drank
2182  use scale_time, only: &
2183  nowsec => time_nowdaysec
2184  implicit none
2185 
2186  integer, intent(in) :: fid
2187  integer, intent(in) :: vid
2188  real(RP), intent(in) :: var(:)
2189  character(len=*), intent(in) :: varname
2190  character(len=*), intent(in) :: axistype
2191 
2192  integer :: dim1_S, dim1_E
2193  integer :: rankidx(2)
2194  integer :: start(1) ! used only when IO_AGGREGATE is .true.
2195  logical :: exec = .true.
2196  !---------------------------------------------------------------------------
2197 
2198  call prof_rapstart('FILE_O_NetCDF', 2)
2199 
2200  rankidx(1) = prc_2drank(prc_myrank,1)
2201  rankidx(2) = prc_2drank(prc_myrank,2)
2202 
2203  if ( axistype == 'Z' ) then
2204  dim1_s = ks
2205  dim1_e = ke
2206  start(1) = 1
2207  if( io_aggregate .AND. prc_myrank > 0 ) exec = .false. ! only rank 0 writes
2208  elseif( axistype == 'X' ) then
2209  dim1_s = isb
2210  dim1_e = ieb
2211  start(1) = isga
2212  if( io_aggregate .AND. rankidx(2) > 0 ) exec = .false. ! only south most row processes write
2213  elseif( axistype == 'Y' ) then
2214  dim1_s = jsb
2215  dim1_e = jeb
2216  start(1) = jsga
2217  if( io_aggregate .AND. rankidx(1) > 0 ) exec = .false. ! only west most column processes write
2218  else
2219  write(*,*) 'xxx [FILEIO_write_var_1D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
2220  call prc_mpistop
2221  endif
2222 
2223  if( exec ) call filewrite( fid, vid, var(dim1_s:dim1_e), &
2224  nowsec, nowsec, start ) ! [IN]
2225 
2226  call prof_rapend ('FILE_O_NetCDF', 2)
2227 
2228  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
integer, public ieb
integer, public jsga
start point of the full domain: cy, global
module TIME
Definition: scale_time.F90:15
module PROCESS
integer, public prc_myrank
process num in local communicator
module RM PROCESS
integer, public isga
start point of the full domain: cx, global
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_write_var_2d()

subroutine scale_fileio::fileio_write_var_2d ( integer, intent(in)  fid,
integer, intent(in)  vid,
real(rp), dimension(:,:), intent(in)  var,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
logical, intent(in), optional  nohalo 
)

Write 2D data to file.

Parameters
[in]fidfile ID
[in]vidnetCDF variable ID
[in]varvalue of the variable
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y)
[in]nohaloswitch whether include halo data or not (default=false)

Definition at line 2240 of file scale_fileio.F90.

References scale_grid_index::ia, scale_grid_index::ie, scale_grid_index::ieb, scale_stdio::io_aggregate, scale_grid_index::is, scale_grid_index::isb, scale_grid_index::isga, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::jeb, scale_grid_index::js, scale_grid_index::jsb, scale_grid_index::jsga, scale_grid_index::ke, scale_grid_index::ks, scale_rm_process::prc_2drank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), gtool_file::rmiss, and scale_time::time_nowdaysec.

Referenced by fileio_write_2d().

2240  use gtool_file, only: &
2241  rmiss
2242  use gtool_file, only: &
2243  filewrite
2244  use scale_process, only: &
2245  prc_myrank, &
2246  prc_mpistop
2247  use scale_rm_process, only: &
2248  prc_2drank, &
2249  prc_num_x, &
2250  prc_num_y
2251  use scale_time, only: &
2252  nowsec => time_nowdaysec
2253  implicit none
2254 
2255  integer, intent(in) :: fid
2256  integer, intent(in) :: vid
2257  real(RP), intent(in) :: var(:,:)
2258  character(len=*), intent(in) :: varname
2259  character(len=*), intent(in) :: axistype
2260  logical, intent(in), optional :: nohalo
2261 
2262  real(RP) :: varhalo( size(var(:,1)), size(var(1,:)) )
2263 
2264  integer :: dim1_S, dim1_E
2265  integer :: dim2_S, dim2_E
2266 
2267  integer :: i, j
2268  logical :: nohalo_
2269  integer :: rankidx(2)
2270  integer :: start(2) ! used only when IO_AGGREGATE is .true.
2271  logical :: exec = .true.
2272  !---------------------------------------------------------------------------
2273 
2274  call prof_rapstart('FILE_O_NetCDF', 2)
2275 
2276  rankidx(1) = prc_2drank(prc_myrank,1)
2277  rankidx(2) = prc_2drank(prc_myrank,2)
2278 
2279  start(1) = isga
2280  start(2) = jsga
2281 
2282  nohalo_ = .false.
2283  if( present(nohalo) ) nohalo_ = nohalo
2284 
2285  if ( axistype == 'XY' &
2286  .OR. axistype == 'UY' &
2287  .OR. axistype == 'XV' &
2288  .OR. axistype == 'UV' ) then
2289  dim1_s = isb
2290  dim1_e = ieb
2291  dim2_s = jsb
2292  dim2_e = jeb
2293  if ( io_aggregate ) then
2294  if( rankidx(1) == 0 ) dim1_s = 1
2295  if( rankidx(1) == prc_num_x - 1 ) dim1_e = ia
2296  if( rankidx(2) == 0 ) dim2_s = 1
2297  if( rankidx(2) == prc_num_y - 1 ) dim2_e = ja
2298  end if
2299  elseif( axistype == 'ZX' ) then
2300  dim1_s = ks
2301  dim1_e = ke
2302  dim2_s = isb
2303  dim2_e = ieb
2304  start(2) = start(1)
2305  start(1) = 1
2306  if ( io_aggregate .AND. rankidx(2) > 0 ) then
2307  exec = .false. ! only south most row processes write
2308  if( rankidx(1) == 0 ) dim2_s = 1
2309  if( rankidx(1) == prc_num_x - 1 ) dim2_e = ia
2310  end if
2311  else
2312  write(*,*) 'xxx [FILEIO_write_var_2D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
2313  call prc_mpistop
2314  endif
2315 
2316  if ( exec ) then
2317  if ( nohalo_ ) then ! fill halo cells with RMISS
2318  varhalo(:,:) = var(:,:)
2319 
2320  ! W halo
2321  do j = 1, ja
2322  do i = 1, is-1
2323  varhalo(i,j) = rmiss
2324  end do
2325  end do
2326  ! E halo
2327  do j = 1, ja
2328  do i = ie+1, ia
2329  varhalo(i,j) = rmiss
2330  end do
2331  end do
2332  ! S halo
2333  do j = 1, js-1
2334  do i = 1, ia
2335  varhalo(i,j) = rmiss
2336  end do
2337  end do
2338  ! N halo
2339  do j = je+1, ja
2340  do i = 1, ia
2341  varhalo(i,j) = rmiss
2342  end do
2343  end do
2344 
2345  call filewrite( fid, vid, varhalo(dim1_s:dim1_e,dim2_s:dim2_e), &
2346  nowsec, nowsec, start ) ! [IN]
2347  else
2348  call filewrite( fid, vid, var(dim1_s:dim1_e,dim2_s:dim2_e), &
2349  nowsec, nowsec, start ) ! [IN]
2350  end if
2351 
2352  end if
2353 
2354  call prof_rapend ('FILE_O_NetCDF', 2)
2355 
2356  return
integer, public prc_num_x
x length of 2D processor topology
integer, public is
start point of inner domain: x, local
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
real(dp), parameter, public rmiss
Definition: gtool_file.f90:150
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
integer, public jsga
start point of the full domain: cy, global
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module PROCESS
integer, public prc_myrank
process num in local communicator
module RM PROCESS
integer, public isga
start point of the full domain: cx, global
integer, public ie
end point of inner domain: x, local
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_write_var_3d()

subroutine scale_fileio::fileio_write_var_3d ( integer, intent(in)  fid,
integer, intent(in)  vid,
real(rp), dimension(:,:,:), intent(in)  var,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
logical, intent(in), optional  nohalo 
)

Write 3D data to file.

Parameters
[in]fidfile ID
[in]vidnetCDF variable ID
[in]varvalue of the variable
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y)
[in]nohaloinclude halo data?

Definition at line 2368 of file scale_fileio.F90.

References scale_grid_index::ia, scale_grid_index::ie, scale_grid_index::ieb, scale_stdio::io_aggregate, scale_grid_index::is, scale_grid_index::isb, scale_grid_index::isga, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::jeb, scale_grid_index::js, scale_grid_index::jsb, scale_grid_index::jsga, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, scale_land_grid_index::lke, scale_land_grid_index::lkmax, scale_land_grid_index::lks, scale_rm_process::prc_2drank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), gtool_file::rmiss, scale_time::time_nowdaysec, scale_urban_grid_index::uke, scale_urban_grid_index::ukmax, and scale_urban_grid_index::uks.

Referenced by fileio_write_3d().

2368  use gtool_file, only: &
2369  rmiss
2370  use gtool_file, only: &
2371  filewrite
2372  use scale_process, only: &
2373  prc_myrank, &
2374  prc_mpistop
2375  use scale_rm_process, only: &
2376  prc_2drank, &
2377  prc_num_x, &
2378  prc_num_y
2379  use scale_time, only: &
2380  nowsec => time_nowdaysec
2381  implicit none
2382 
2383  integer, intent(in) :: fid
2384  integer, intent(in) :: vid
2385  real(RP), intent(in) :: var(:,:,:)
2386  character(len=*), intent(in) :: varname
2387  character(len=*), intent(in) :: axistype
2388 
2389  logical, intent(in), optional :: nohalo
2390 
2391  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)), size(var(1,1,:)) )
2392 
2393  integer :: dim1_S, dim1_E, dim1_max
2394  integer :: dim2_S, dim2_E
2395  integer :: dim3_S, dim3_E
2396 
2397  integer :: i, j, k
2398  logical :: nohalo_
2399  integer :: rankidx(2)
2400  integer :: start(3) ! used only when IO_AGGREGATE is .true.
2401  !---------------------------------------------------------------------------
2402 
2403  call prof_rapstart('FILE_O_NetCDF', 2)
2404 
2405  nohalo_ = .false.
2406  if( present(nohalo) ) nohalo_ = nohalo
2407 
2408  rankidx(1) = prc_2drank(prc_myrank,1)
2409  rankidx(2) = prc_2drank(prc_myrank,2)
2410 
2411  start(1) = 1
2412  start(2) = isga
2413  start(3) = jsga
2414 
2415  dim2_s = isb
2416  dim2_e = ieb
2417  dim3_s = jsb
2418  dim3_e = jeb
2419  if ( io_aggregate ) then
2420  if( rankidx(1) == 0 ) dim2_s = 1
2421  if( rankidx(1) == prc_num_x - 1 ) dim2_e = ia
2422  if( rankidx(2) == 0 ) dim3_s = 1
2423  if( rankidx(2) == prc_num_y - 1 ) dim3_e = ja
2424  end if
2425 
2426  if ( axistype == 'ZXY' &
2427  .OR. axistype == 'ZHXY' &
2428  .OR. axistype == 'ZXHY' &
2429  .OR. axistype == 'ZXYH' ) then
2430  dim1_max = kmax
2431  dim1_s = ks
2432  dim1_e = ke
2433  elseif( axistype == 'Land' ) then
2434  dim1_max = lkmax
2435  dim1_s = lks
2436  dim1_e = lke
2437  elseif( axistype == 'Urban' ) then
2438  dim1_max = ukmax
2439  dim1_s = uks
2440  dim1_e = uke
2441  else
2442  write(*,*) 'xxx [FILEIO_write_var_3D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
2443  call prc_mpistop
2444  endif
2445 
2446  if ( nohalo_ ) then
2447  varhalo(:,:,:) = var(:,:,:)
2448 
2449  ! W halo
2450  do k = 1, dim1_max
2451  do j = 1, ja
2452  do i = 1, is-1
2453  varhalo(k,i,j) = rmiss
2454  end do
2455  end do
2456  end do
2457  ! E halo
2458  do k = 1, dim1_max
2459  do j = 1, ja
2460  do i = ie+1, ia
2461  varhalo(k,i,j) = rmiss
2462  end do
2463  end do
2464  end do
2465  ! S halo
2466  do k = 1, dim1_max
2467  do j = 1, js-1
2468  do i = 1, ia
2469  varhalo(k,i,j) = rmiss
2470  end do
2471  end do
2472  end do
2473  ! N halo
2474  do k = 1, dim1_max
2475  do j = je+1, ja
2476  do i = 1, ia
2477  varhalo(k,i,j) = rmiss
2478  end do
2479  end do
2480  end do
2481 
2482  call filewrite( fid, vid, varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e), &
2483  nowsec, nowsec, start ) ! [IN]
2484  else
2485  call filewrite( fid, vid, var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e), &
2486  nowsec, nowsec, start ) ! [IN]
2487  end if
2488 
2489  call prof_rapend ('FILE_O_NetCDF', 2)
2490 
2491  return
integer, public prc_num_x
x length of 2D processor topology
integer, public is
start point of inner domain: x, local
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
real(dp), parameter, public rmiss
Definition: gtool_file.f90:150
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
integer, public jsga
start point of the full domain: cy, global
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module PROCESS
integer, public prc_myrank
process num in local communicator
module RM PROCESS
integer, public isga
start point of the full domain: cx, global
integer, public ie
end point of inner domain: x, local
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_write_var_3d_t()

subroutine scale_fileio::fileio_write_var_3d_t ( integer, intent(in)  fid,
integer, intent(in)  vid,
real(rp), dimension(:,:,:), intent(in)  var,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
real(rp), intent(in)  timeintv,
integer, intent(in), optional  timetarg,
logical, intent(in), optional  nohalo 
)

Write 3D data with time dimension to file.

Parameters
[in]fidfile ID
[in]vidnetCDF variable ID
[in]varvalue of the variable
[in]varnamename of the variable
[in]axistypeaxis type (X/Y/Time)
[in]timeintvtime interval [sec]
[in]timetargtarget timestep (optional)
[in]nohaloinclude halo data?

Definition at line 2505 of file scale_fileio.F90.

References scale_grid_index::ia, scale_grid_index::ie, scale_grid_index::ieb, scale_stdio::io_aggregate, scale_grid_index::is, scale_grid_index::isb, scale_grid_index::isga, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::jeb, scale_grid_index::js, scale_grid_index::jsb, scale_grid_index::jsga, scale_rm_process::prc_2drank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and gtool_file::rmiss.

Referenced by fileio_write_3d_t().

2505  use gtool_file, only: &
2506  rmiss
2507  use gtool_file, only: &
2508  filewrite
2509  use scale_process, only: &
2510  prc_myrank, &
2511  prc_mpistop
2512  use scale_rm_process, only: &
2513  prc_2drank, &
2514  prc_num_x, &
2515  prc_num_y
2516  implicit none
2517 
2518  integer, intent(in) :: fid
2519  integer, intent(in) :: vid
2520  real(RP), intent(in) :: var(:,:,:)
2521  character(len=*), intent(in) :: varname
2522  character(len=*), intent(in) :: axistype
2523  real(RP), intent(in) :: timeintv
2524 
2525  integer, intent(in), optional :: timetarg
2526  logical, intent(in), optional :: nohalo
2527 
2528  real(RP) :: varhalo( size(var(:,1,1)), size(var(1,:,1)) )
2529 
2530  integer :: dim1_S, dim1_E
2531  integer :: dim2_S, dim2_E
2532 
2533  real(DP) :: time_interval, nowtime
2534 
2535  integer :: step
2536  integer :: i, j, n
2537  logical :: nohalo_
2538  integer :: rankidx(2)
2539  integer :: start(3) ! used only when IO_AGGREGATE is .true.
2540  !---------------------------------------------------------------------------
2541 
2542  call prof_rapstart('FILE_O_NetCDF', 2)
2543 
2544  nohalo_ = .false.
2545  if( present(nohalo) ) nohalo_ = nohalo
2546 
2547  time_interval = timeintv
2548  step = size(var(isb,jsb,:))
2549 
2550  rankidx(1) = prc_2drank(prc_myrank,1)
2551  rankidx(2) = prc_2drank(prc_myrank,2)
2552 
2553  if ( axistype == 'XYT' ) then
2554  dim1_s = isb
2555  dim1_e = ieb
2556  dim2_s = jsb
2557  dim2_e = jeb
2558  if ( io_aggregate ) then
2559  if( rankidx(1) == 0 ) dim1_s = 1
2560  if( rankidx(1) == prc_num_x - 1 ) dim1_e = ia
2561  if( rankidx(2) == 0 ) dim2_s = 1
2562  if( rankidx(2) == prc_num_y - 1 ) dim2_e = ja
2563  end if
2564  else
2565  write(*,*) 'xxx [FILEIO_write_var_3D_t] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
2566  call prc_mpistop
2567  endif
2568 
2569  start(1) = isga
2570  start(2) = jsga
2571  ! start(3) time dimension will be set in file_write_data()
2572 
2573  if ( present(timetarg) ) then
2574  nowtime = (timetarg-1) * time_interval
2575 
2576  if ( nohalo_ ) then
2577  varhalo(:,:) = var(:,:,timetarg)
2578 
2579  ! W halo
2580  do j = 1, ja
2581  do i = 1, is-1
2582  varhalo(i,j) = rmiss
2583  end do
2584  end do
2585  ! E halo
2586  do j = 1, ja
2587  do i = ie+1, ia
2588  varhalo(i,j) = rmiss
2589  end do
2590  end do
2591  ! S halo
2592  do j = 1, js-1
2593  do i = 1, ia
2594  varhalo(i,j) = rmiss
2595  end do
2596  end do
2597  ! N halo
2598  do j = je+1, ja
2599  do i = 1, ia
2600  varhalo(i,j) = rmiss
2601  end do
2602  end do
2603 
2604  call filewrite( fid, vid, varhalo(dim1_s:dim1_e,dim2_s:dim2_e), &
2605  nowtime, nowtime, start ) ! [IN]
2606  else
2607  call filewrite( fid, vid, var(dim1_s:dim1_e,dim2_s:dim2_e,timetarg), &
2608  nowtime, nowtime, start ) ! [IN]
2609  end if
2610  else
2611  nowtime = 0.0_dp
2612  do n = 1, step
2613  if ( nohalo_ ) then
2614  varhalo(:,:) = var(:,:,n)
2615 
2616  ! W halo
2617  do j = 1, ja
2618  do i = 1, is-1
2619  varhalo(i,j) = rmiss
2620  end do
2621  end do
2622  ! E halo
2623  do j = 1, ja
2624  do i = ie+1, ia
2625  varhalo(i,j) = rmiss
2626  end do
2627  end do
2628  ! S halo
2629  do j = 1, js-1
2630  do i = 1, ia
2631  varhalo(i,j) = rmiss
2632  end do
2633  end do
2634  ! N halo
2635  do j = je+1, ja
2636  do i = 1, ia
2637  varhalo(i,j) = rmiss
2638  end do
2639  end do
2640 
2641  call filewrite( fid, vid, varhalo(dim1_s:dim1_e,dim2_s:dim2_e), &
2642  nowtime, nowtime, start ) ! [IN]
2643  else
2644  call filewrite( fid, vid, var(dim1_s:dim1_e,dim2_s:dim2_e,n), &
2645  nowtime, nowtime, start ) ! [IN]
2646  end if
2647  nowtime = nowtime + time_interval
2648  enddo
2649  endif
2650 
2651  call prof_rapend ('FILE_O_NetCDF', 2)
2652 
2653  return
integer, public prc_num_x
x length of 2D processor topology
integer, public is
start point of inner domain: x, local
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
real(dp), parameter, public rmiss
Definition: gtool_file.f90:150
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
integer, public jsga
start point of the full domain: cy, global
integer, public js
start point of inner domain: y, local
module PROCESS
integer, public prc_myrank
process num in local communicator
module RM PROCESS
integer, public isga
start point of the full domain: cx, global
integer, public ie
end point of inner domain: x, local
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ fileio_write_var_4d()

subroutine scale_fileio::fileio_write_var_4d ( integer, intent(in)  fid,
integer, intent(in)  vid,
real(rp), dimension(:,:,:,:), intent(in)  var,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  axistype,
real(rp), intent(in)  timeintv,
integer, intent(in), optional  timetarg,
logical, intent(in), optional  nohalo 
)

Write 4D data to file.

Parameters
[in]fidfile ID
[in]vidnetCDF variable ID
[in]varvalue of the variable
[in]varnamename of the variable
[in]axistypeaxis type (Z/X/Y/Time)
[in]timeintvtime interval [sec]
[in]timetargtarget timestep (optional)
[in]nohaloinclude halo data?

Definition at line 2667 of file scale_fileio.F90.

References scale_const::const_eps, fileio_close(), scale_grid_index::ia, scale_grid_index::ie, scale_grid_index::ieb, scale_stdio::io_aggregate, scale_grid_index::is, scale_grid_index::isb, scale_grid_index::isga, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::jeb, scale_grid_index::js, scale_grid_index::jsb, scale_grid_index::jsga, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, scale_rm_process::prc_2drank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and gtool_file::rmiss.

Referenced by fileio_write_4d().

2667  use gtool_file, only: &
2668  rmiss
2669  use gtool_file, only: &
2670  filewrite
2671  use scale_process, only: &
2672  prc_myrank, &
2673  prc_mpistop
2674  use scale_rm_process, only: &
2675  prc_2drank, &
2676  prc_num_x, &
2677  prc_num_y
2678  implicit none
2679 
2680  integer, intent(in) :: fid
2681  integer, intent(in) :: vid
2682  real(RP), intent(in) :: var(:,:,:,:)
2683  character(len=*), intent(in) :: varname
2684  character(len=*), intent(in) :: axistype
2685  real(RP), intent(in) :: timeintv
2686 
2687  integer, intent(in), optional :: timetarg
2688  logical, intent(in), optional :: nohalo
2689 
2690  real(RP) :: varhalo( size(var(:,1,1,1)), size(var(1,:,1,1)), size(var(1,1,:,1)) )
2691 
2692  integer :: dim1_S, dim1_E, dim1_max
2693  integer :: dim2_S, dim2_E
2694  integer :: dim3_S, dim3_E
2695 
2696  real(DP) :: time_interval, nowtime
2697 
2698  integer :: step
2699  integer :: i, j, k, n
2700  logical :: nohalo_
2701  integer :: rankidx(2)
2702  integer :: start(4) ! used only when IO_AGGREGATE is .true.
2703  !---------------------------------------------------------------------------
2704 
2705  call prof_rapstart('FILE_O_NetCDF', 2)
2706 
2707  nohalo_ = .false.
2708  if( present(nohalo) ) nohalo_ = nohalo
2709 
2710  rankidx(1) = prc_2drank(prc_myrank,1)
2711  rankidx(2) = prc_2drank(prc_myrank,2)
2712 
2713  start(1) = 1
2714  start(2) = isga
2715  start(3) = jsga
2716  ! start(4) time dimension will be set in file_write_data()
2717 
2718  time_interval = timeintv
2719  step = size(var,4)
2720 
2721  if ( axistype == 'ZXYT' ) then
2722  dim1_max = kmax
2723  dim1_s = ks
2724  dim1_e = ke
2725  dim2_s = isb
2726  dim2_e = ieb
2727  dim3_s = jsb
2728  dim3_e = jeb
2729  if ( io_aggregate ) then
2730  if( rankidx(1) == 0 ) dim2_s = 1
2731  if( rankidx(1) == prc_num_x - 1 ) dim2_e = ia
2732  if( rankidx(2) == 0 ) dim3_s = 1
2733  if( rankidx(2) == prc_num_y - 1 ) dim3_e = ja
2734  end if
2735  else
2736  write(*,*) 'xxx [FILEIO_write_var_4D] unsupported axis type. Check! axistype:', trim(axistype), ', item:',trim(varname)
2737  call prc_mpistop
2738  endif
2739 
2740  if ( present(timetarg) ) then
2741  nowtime = (timetarg-1) * time_interval
2742 
2743  if ( nohalo_ ) then
2744  varhalo(:,:,:) = var(:,:,:,timetarg)
2745 
2746  ! W halo
2747  do k = 1, dim1_max
2748  do j = 1, ja
2749  do i = 1, is-1
2750  varhalo(k,i,j) = rmiss
2751  end do
2752  end do
2753  end do
2754  ! E halo
2755  do k = 1, dim1_max
2756  do j = 1, ja
2757  do i = ie+1, ia
2758  varhalo(k,i,j) = rmiss
2759  end do
2760  end do
2761  end do
2762  ! S halo
2763  do k = 1, dim1_max
2764  do j = 1, js-1
2765  do i = 1, ia
2766  varhalo(k,i,j) = rmiss
2767  end do
2768  end do
2769  end do
2770  ! N halo
2771  do k = 1, dim1_max
2772  do j = je+1, ja
2773  do i = 1, ia
2774  varhalo(k,i,j) = rmiss
2775  end do
2776  end do
2777  end do
2778 
2779  call filewrite( fid, vid, varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e), &
2780  nowtime, nowtime, start ) ! [IN]
2781  else
2782  call filewrite( fid, vid, var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e,timetarg), &
2783  nowtime, nowtime, start ) ! [IN]
2784  end if
2785  else
2786  nowtime = 0.0_dp
2787  do n = 1, step
2788  if ( nohalo_ ) then
2789  varhalo(:,:,:) = var(:,:,:,n)
2790 
2791  ! W halo
2792  do k = 1, dim1_max
2793  do j = 1, ja
2794  do i = 1, is-1
2795  varhalo(k,i,j) = rmiss
2796  end do
2797  end do
2798  end do
2799  ! E halo
2800  do k = 1, dim1_max
2801  do j = 1, ja
2802  do i = ie+1, ia
2803  varhalo(k,i,j) = rmiss
2804  end do
2805  end do
2806  end do
2807  ! S halo
2808  do k = 1, dim1_max
2809  do j = 1, js-1
2810  do i = 1, ia
2811  varhalo(k,i,j) = rmiss
2812  end do
2813  end do
2814  end do
2815  ! N halo
2816  do k = 1, dim1_max
2817  do j = je+1, ja
2818  do i = 1, ia
2819  varhalo(k,i,j) = rmiss
2820  end do
2821  end do
2822  end do
2823 
2824  call filewrite( fid, vid, varhalo(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e), &
2825  nowtime, nowtime, start ) ! [IN]
2826  else
2827  call filewrite( fid, vid, var(dim1_s:dim1_e,dim2_s:dim2_e,dim3_s:dim3_e,n), &
2828  nowtime, nowtime, start ) ! [IN]
2829  end if
2830  nowtime = nowtime + time_interval
2831  enddo
2832  endif
2833 
2834  call prof_rapend ('FILE_O_NetCDF', 2)
2835 
2836  return
integer, public prc_num_x
x length of 2D processor topology
integer, public is
start point of inner domain: x, local
module GTOOL_FILE
Definition: gtool_file.f90:17
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
real(dp), parameter, public rmiss
Definition: gtool_file.f90:150
integer, public ieb
integer, public prc_num_y
y length of 2D processor topology
integer, public jsga
start point of the full domain: cy, global
integer, public js
start point of inner domain: y, local
module PROCESS
integer, public prc_myrank
process num in local communicator
module RM PROCESS
integer, public isga
start point of the full domain: cx, global
integer, public ie
end point of inner domain: x, local
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, public isb
integer, public jsb
Here is the call graph for this function:
Here is the caller graph for this function: