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

module file_grads More...

Functions/Subroutines

subroutine, public file_grads_open (file_name, file_id)
 Open. More...
 
subroutine, public file_grads_varid (file_id, var_name, var_id)
 
logical function, public file_grads_isoned (file_id, var_id)
 
subroutine file_grads_get_shape_name (file_id, var_name, shape)
 
subroutine file_grads_get_shape_id (file_id, var_id, shape)
 
subroutine file_grads_read_1d_id (file_id, var_id, var, step, postfix)
 
subroutine file_grads_read_2d_id (file_id, var_id, var, step, postfix)
 
subroutine file_grads_read_3d_id (file_id, var_id, var, step, postfix)
 
subroutine, public file_grads_close (file_id)
 
subroutine read_data_int1 (fid, irecl, nx, ny, nz, k, yrev, var, ierr)
 
character(len=len(str)) function upcase (str)
 

Detailed Description

module file_grads

Description
read data from GrADS file
Author
Team SCALE
NAMELIST
  • GrADS_DIMS
    nametypedefault valuecomment
    NX integer
    NY integer
    NZ integer

  • GrADS_ITEM
    nametypedefault valuecomment
    NAME character(len=H_SHORT) up to 16 characters
    DTYPE character(len=H_SHORT) 'linear','levels','map'
    FNAME character(len=H_LONG) head of file name
    SWPOINT real(RP) start point (south-west point) for "linear"
    DD real(RP) dlon,dlat for "linear"
    LNUM integer number of data
    LVARS real(RP), dimension(LVARS_MAX) values for "levels"
    STARTREC integer record position
    TOTALREC integer total record number per one time
    MISSVAL real(SP) missing value
    NX integer
    NY integer
    NZ integer
    YREV logical
    BINTYPE character(len=H_SHORT) binary type: 'int?' or 'real?' ?=2,4, or 8

  • nml_grads_grid
    nametypedefault valuecomment
    DUMMY logical

  • grdvar
    nametypedefault valuecomment
    DUMMY logical

History Output
No history output

Function/Subroutine Documentation

◆ file_grads_open()

subroutine, public scale_file_grads::file_grads_open ( character(len=*), intent(in)  file_name,
integer, intent(out)  file_id 
)

Open.

Definition at line 102 of file scale_file_grads.F90.

102  use scale_prc, only: &
103  prc_abort
104  use scale_const, only: &
105  undef => const_undef
106  implicit none
107  character(len=*), intent(in) :: file_name
108  integer, intent(out) :: file_id
109 
110  character(len=H_SHORT) :: name ! up to 16 characters
111  character(len=H_SHORT) :: dtype ! 'linear','levels','map'
112  character(len=H_LONG) :: fname ! head of file name
113  real(RP) :: swpoint ! start point (south-west point) for "linear"
114  real(RP) :: dd ! dlon,dlat for "linear"
115  integer :: lnum ! number of data
116  real(RP) :: lvars(lvars_max) ! values for "levels"
117  integer :: startrec ! record position
118  integer :: totalrec ! total record number per one time
119  real(SP) :: missval ! missing value
120  integer :: nx ! optional
121  integer :: ny ! optional
122  integer :: nz ! optional
123  character(len=H_SHORT) :: fendian ! option for "map"
124  logical :: yrev ! option for "map", if yrev=.true., order of data is NW to SE.
125  character(len=H_SHORT) :: bintype ! binary type: 'int?' or 'real?' ?=2,4, or 8
126 
127  namelist /grads_dims/ &
128  nx, &
129  ny, &
130  nz
131 
132  namelist /grads_item/ &
133  name, & ! necessary
134  dtype, & ! necessary
135  fname, & ! necessary except for linear data
136  swpoint, & ! for linear data
137  dd, & ! for linear data
138  lnum, & ! for levels data
139  lvars, & ! for levels data
140  startrec, & ! for map data
141  totalrec, & ! for map data
142  missval, & ! option
143  nx, & ! option
144  ny, & ! option
145  nz, & ! option
146  yrev, & ! option
147  bintype ! option
148 ! fendian ! option
149 
150  character(len=H_LONG) :: dirname
151 
152  integer :: fid
153  integer :: nvars
154  integer :: ierr
155  integer :: n
156  !---------------------------------------------------------------------------
157 
158  log_newline
159  log_info("FILE_GrADS_open",*) 'open namelist file :', trim(file_name)
160 
161  ! check exist
162  do n = 1, nnmls
163  if ( nmls(n)%fname == file_name ) then
164  ! alread read
165  file_id = n
166  return
167  end if
168  end do
169 
170 
171  fid = io_get_available_fid()
172  !--- open namelist
173  open( fid, &
174  file = trim(file_name), &
175  form = 'formatted', &
176  status = 'old', &
177  action = 'read', &
178  iostat = ierr )
179  if ( ierr /= 0 ) then
180  log_error("FILE_GrADS_open",*) 'Input file is not found! ', trim(file_name)
181  call prc_abort
182  endif
183 
184  call check_oldnamelist( fid )
185 
186  !--- read namelist dims
187  read(fid,nml=grads_dims,iostat=ierr)
188  if( ierr /= 0 ) then !--- missing or fatal error
189  log_error("FILE_GrADS_open",*) 'Not appropriate names in GrADS_DIMS in ', trim(file_name),'. Check!'
190  call prc_abort
191  endif
192  log_nml(grads_dims)
193 
194 
195  nnmls = nnmls + 1
196  file_id = nnmls
197 
198  nmls(file_id)%fname = file_name
199  nmls(file_id)%nx = nx
200  nmls(file_id)%ny = ny
201  nmls(file_id)%nz = nz
202 
203 
204 
205  !--- count the number of variables
206  rewind(fid)
207  nvars = 0
208  do n = 1, vars_max
209  read(fid, nml=grads_item, iostat=ierr)
210  if( ierr > 0 )then
211  log_error("FILE_GrADS_open",*) 'Not appropriate names in GrADS_ITEM in ', trim(file_name),'. Check!'
212  call prc_abort
213  else if( ierr < 0 )then
214  exit
215  endif
216  nvars = nvars + 1
217  enddo
218 
219  if ( nvars > vars_max ) then
220  log_error("FILE_GRADS_open",*) 'The number of grads vars exceeds the limit! ', &
221  nvars, ' > ', vars_max
222  call prc_abort
223  endif
224 
225  nmls(file_id)%nvars = nvars
226  allocate( nmls(file_id)%vars(nvars) )
227 
228  n = index( file_name, '/', back=.true. )
229  if ( n > 0 ) then
230  dirname = file_name(1:n)
231  else
232  dirname = ""
233  end if
234 
235  !--- read information of the variables
236  rewind(fid)
237  do n = 1, nvars
238 
239  ! set default
240  name = ''
241  dtype = ''
242  fname = ''
243  swpoint = undef
244  dd = undef
245  lnum = -1
246  lvars(:) = undef
247  startrec = -1
248  totalrec = -1
249  nx = nmls(file_id)%nx
250  ny = nmls(file_id)%ny
251  nz = nmls(file_id)%nz
252  yrev = .false.
253  fendian = 'big'
254  missval = undef
255  bintype = 'real4'
256 
257  ! read namelist
258  read(fid, nml=grads_item, iostat=ierr)
259  if( ierr /= 0 ) exit
260 
261  nmls(file_id)%vars(n)%name = upcase(name)
262  if ( fname(1:1) == "/" ) then
263  nmls(file_id)%vars(n)%fname = fname
264  else
265  nmls(file_id)%vars(n)%fname = trim(dirname) // fname
266  end if
267  nmls(file_id)%vars(n)%dtype = dtype
268  nmls(file_id)%vars(n)%swpoint = swpoint
269  nmls(file_id)%vars(n)%dd = dd
270  nmls(file_id)%vars(n)%lnum = lnum
271  if ( lnum > 0 ) then
272  allocate( nmls(file_id)%vars(n)%lvars(lnum) )
273  nmls(file_id)%vars(n)%lvars(:) = lvars(1:lnum)
274  end if
275  nmls(file_id)%vars(n)%startrec = startrec
276  nmls(file_id)%vars(n)%totalrec = totalrec
277  nmls(file_id)%vars(n)%missval = missval
278  nmls(file_id)%vars(n)%nx = nx
279  nmls(file_id)%vars(n)%ny = ny
280  nmls(file_id)%vars(n)%nz = nz
281  nmls(file_id)%vars(n)%yrev = yrev
282  if ( fendian == "big" ) then
283  nmls(file_id)%vars(n)%endian = 1
284  else
285  nmls(file_id)%vars(n)%endian = 0
286  end if
287  nmls(file_id)%vars(n)%bintype = bintype
288 
289  end do
290 
291  close( fid )
292 

References scale_const::const_undef, scale_io::io_get_available_fid(), scale_prc::prc_abort(), and upcase().

Referenced by mod_cnv2d::cnv2d_grads_init(), mod_cnvuser::cnvuser(), mod_cnvuser::cnvuser_setup(), mod_copytopo::copytopo_get_data_grads(), mod_copytopo::copytopo_get_size_grads(), mod_copytopo::copytopo_get_size_wrfarw(), mod_realinput_grads::parentatmossetupgrads(), mod_realinput_grads::parentlandsetupgrads(), and mod_realinput_grads::parentoceansetupgrads().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_grads_varid()

subroutine, public scale_file_grads::file_grads_varid ( integer, intent(in)  file_id,
character(len=*), intent(in)  var_name,
integer, intent(out)  var_id 
)

Definition at line 300 of file scale_file_grads.F90.

300  use scale_prc, only: &
301  prc_abort
302  implicit none
303  integer, intent(in) :: file_id
304  character(len=*), intent(in) :: var_name
305  integer, intent(out) :: var_id
306 
307  character(len=len(var_name)) :: vname
308  integer :: n
309 
310  if ( file_id < 0 ) then
311  log_error("FILE_GrADS_varid",*) 'file_id is invalid: ', file_id
312  call prc_abort
313  end if
314 
315  vname = upcase(var_name)
316 
317  var_id = -1
318  do n = 1, nmls(file_id)%nvars
319  if ( nmls(file_id)%vars(n)%name == vname ) then
320  var_id = n
321  return
322  end if
323  end do
324 
325  return

References scale_prc::prc_abort(), and upcase().

Referenced by mod_cnv2d::cnv2d_grads_init(), mod_cnvuser::cnvuser(), mod_cnvuser::cnvuser_setup(), mod_copytopo::copytopo_get_data_grads(), file_grads_get_shape_id(), file_grads_get_shape_name(), file_grads_read_1d_id(), file_grads_read_2d_id(), mod_realinput_grads::parentatmossetupgrads(), mod_realinput_grads::parentlandsetupgrads(), and mod_realinput_grads::parentoceansetupgrads().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_grads_isoned()

logical function, public scale_file_grads::file_grads_isoned ( integer, intent(in)  file_id,
integer, intent(in)  var_id 
)

Definition at line 332 of file scale_file_grads.F90.

332  use scale_prc, only: &
333  prc_abort
334  implicit none
335  integer, intent(in) :: file_id
336  integer, intent(in) :: var_id
337  logical :: ret
338 
339  if ( file_id < 0 ) then
340  log_error("FILE_GrADS_isOneD",*) 'file_id is invalid: ', file_id
341  call prc_abort
342  end if
343  if ( var_id < 0 ) then
344  log_error("FILE_GrADS_isOneD",*) 'var_id is invalid: ', var_id
345  call prc_abort
346  end if
347 
348  select case( nmls(file_id)%vars(var_id)%dtype )
349  case ('linear', 'levels')
350  ret = .true.
351  case default
352  ret = .false.
353  end select
354 
355  return

References scale_prc::prc_abort().

Referenced by mod_cnv2d::cnv2d_grads_init(), mod_cnvuser::cnvuser_setup(), mod_copytopo::copytopo_get_data_grads(), file_grads_get_shape_id(), mod_realinput_grads::parentatmosinputgrads(), mod_realinput_grads::parentatmosopengrads(), mod_realinput_grads::parentlandinputgrads(), and mod_realinput_grads::parentoceaninputgrads().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_grads_get_shape_name()

subroutine scale_file_grads::file_grads_get_shape_name ( integer, intent(in)  file_id,
character(len=*), intent(in)  var_name,
integer, dimension(:), intent(out)  shape 
)

Definition at line 363 of file scale_file_grads.F90.

363  use scale_prc, only: &
364  prc_abort
365  implicit none
366  integer, intent(in) :: file_id
367  character(len=*), intent(in) :: var_name
368  integer, intent(out) :: shape(:)
369 
370  integer :: var_id
371 
372  call file_grads_varid( file_id, var_name, & ! (in)
373  var_id ) ! (out)
374 
375  if ( var_id < 0 ) then
376  log_error("FILE_GrADS_get_shape",*) 'variable "', trim(var_name), ' is not founed in file "', trim(nmls(file_id)%fname), '"'
377  call prc_abort
378  end if
379 
380  call file_grads_get_shape_id( file_id, var_id, & ! (in)
381  shape(:) ) ! (out)
382 
383  return

References file_grads_get_shape_id(), file_grads_varid(), and scale_prc::prc_abort().

Here is the call graph for this function:

◆ file_grads_get_shape_id()

subroutine scale_file_grads::file_grads_get_shape_id ( integer, intent(in)  file_id,
integer, intent(in)  var_id,
integer, dimension(:), intent(out)  shape 
)

Definition at line 390 of file scale_file_grads.F90.

390  implicit none
391  integer, intent(in) :: file_id
392  integer, intent(in) :: var_id
393  integer, intent(out) :: shape(:)
394 
395  intrinsic :: size
396 
397  if ( file_grads_isoned( file_id, var_id ) ) then
398  if ( nmls(file_id)%vars(var_id)%dtype == "levels" ) then
399  shape(1) = nmls(file_id)%vars(var_id)%lnum
400  else
401  shape(1) = -1
402  end if
403  else if ( size(shape) == 2 ) then
404  shape(1) = nmls(file_id)%vars(var_id)%nx
405  shape(2) = nmls(file_id)%vars(var_id)%ny
406  else
407  shape(1) = nmls(file_id)%vars(var_id)%nz
408  shape(2) = nmls(file_id)%vars(var_id)%nx
409  shape(3) = nmls(file_id)%vars(var_id)%ny
410  end if
411 
412  return

References file_grads_isoned(), file_grads_read_1d_id(), file_grads_varid(), and scale_prc::prc_abort().

Referenced by file_grads_get_shape_name().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_grads_read_1d_id()

subroutine scale_file_grads::file_grads_read_1d_id ( integer, intent(in)  file_id,
integer, intent(in)  var_id,
real(rp), dimension(:), intent(out)  var,
integer, intent(in), optional  step,
character(len=*), intent(in), optional  postfix 
)

Definition at line 457 of file scale_file_grads.F90.

457  implicit none
458  integer, intent(in) :: file_id
459  integer, intent(in) :: var_id
460  real(RP), intent(out) :: var(:)
461  integer, intent(in), optional :: step
462  character(len=*), intent(in), optional :: postfix
463 
464  logical :: exist
465  integer :: vid
466 
467  intrinsic :: size
468  !---------------------------------------------------------------------------
469 
470  if ( file_id < 0 ) then
471  log_error("FILE_GrADS_read_1D_vid",*) 'file_id is invalid: ', file_id
472  end if
473  if ( var_id < 0 ) then
474  log_error("FILE_GrADS_read_1D_vid",*) 'var_id is invalid: ', var_id
475  end if
476 
477  call file_grads_read_data( nmls(file_id)%vars(var_id), & ! (in)
478  1, size(var), & ! (in)
479  var(:), & ! (out)
480  step = step, & ! (int)
481  postfix = postfix ) ! (in)
482 
483  return

References file_grads_read_2d_id(), file_grads_varid(), and scale_prc::prc_abort().

Referenced by file_grads_get_shape_id().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_grads_read_2d_id()

subroutine scale_file_grads::file_grads_read_2d_id ( integer, intent(in)  file_id,
integer, intent(in)  var_id,
real(rp), dimension(:,:), intent(out)  var,
integer, intent(in), optional  step,
character(len=*), intent(in), optional  postfix 
)

Definition at line 527 of file scale_file_grads.F90.

527  implicit none
528  integer, intent(in) :: file_id
529  integer, intent(in) :: var_id
530  real(RP), intent(out) :: var(:,:)
531  integer, intent(in), optional :: step
532  character(len=*), intent(in), optional :: postfix
533 
534  integer :: vid
535  intrinsic :: size
536  !---------------------------------------------------------------------------
537 
538  if ( file_id < 0 ) then
539  log_error("FILE_GrADS_read_2D_vid",*) 'file_id is invalid: ', file_id
540  end if
541  if ( var_id < 0 ) then
542  log_error("FILE_GrADS_read_2D_vid",*) 'var_id is invalid: ', var_id
543  end if
544 
545  call file_grads_read_data( nmls(file_id)%vars(var_id), & ! (in)
546  2, size(var), & ! (in)
547  var(:,:), & ! (out)
548  step = step, & ! (in)
549  postfix = postfix ) ! (in)
550 
551  return

References file_grads_read_3d_id(), file_grads_varid(), and scale_prc::prc_abort().

Referenced by file_grads_read_1d_id().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_grads_read_3d_id()

subroutine scale_file_grads::file_grads_read_3d_id ( integer, intent(in)  file_id,
integer, intent(in)  var_id,
real(rp), dimension(:,:,:), intent(out)  var,
integer, intent(in), optional  step,
character(len=*), intent(in), optional  postfix 
)

Definition at line 595 of file scale_file_grads.F90.

595  implicit none
596  integer, intent(in) :: file_id
597  integer, intent(in) :: var_id
598  real(RP), intent(out) :: var(:,:,:)
599  integer, intent(in), optional :: step
600  character(len=*), intent(in), optional :: postfix
601 
602  integer :: vid
603  intrinsic :: size
604  !---------------------------------------------------------------------------
605 
606  if ( file_id < 0 ) then
607  log_error("FILE_GrADS_read_3D_vid",*) 'file_id is invalid: ', file_id
608  end if
609  if ( var_id < 0 ) then
610  log_error("FILE_GrADS_read_3D_vid",*) 'var_id is invalid: ', var_id
611  end if
612 
613  call file_grads_read_data( nmls(file_id)%vars(var_id), & ! (in)
614  3, size(var), & ! (in)
615  var(:,:,:), & ! (out)
616  step = step, postfix = postfix ) ! (in)
617 
618  return

Referenced by file_grads_read_2d_id().

Here is the caller graph for this function:

◆ file_grads_close()

subroutine, public scale_file_grads::file_grads_close ( integer, intent(in)  file_id)

Definition at line 624 of file scale_file_grads.F90.

624  implicit none
625  integer, intent(in) :: file_id
626 
627  integer :: n, m
628 
629  if ( file_id < 0 ) return
630 
631  do n = 1, nmls(file_id)%nvars
632  do m = 1, nfiles
633  if ( files(m)%fname == nmls(file_id)%vars(n)%fname ) then
634  if ( files(m)%fid > 0 ) then
635  close( files(m)%fid )
636  files(m)%fid = -1
637  files(m)%postfix = ""
638  end if
639  exit
640  end if
641  end do
642  if ( nmls(file_id)%vars(n)%lnum > 0 ) deallocate( nmls(file_id)%vars(n)%lvars )
643  nmls(file_id)%vars(n)%lnum = -1
644  end do
645  deallocate( nmls(file_id)%vars )
646  nmls(file_id)%fname = ""
647  nmls(file_id)%nvars = 0
648 
649  return

References scale_const::const_eps, scale_const::const_undef, scale_io::io_get_available_fid(), scale_prc::prc_abort(), and read_data_int1().

Referenced by mod_copytopo::copytopo_get_data_grads().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_data_int1()

subroutine scale_file_grads::read_data_int1 ( integer, intent(in)  fid,
integer, intent(in)  irecl,
integer, intent(in)  nx,
integer, intent(in)  ny,
integer, intent(in)  nz,
integer, intent(in)  k,
logical, intent(in)  yrev,
real(rp), dimension(:), intent(out)  var,
integer, intent(out)  ierr 
)

Definition at line 863 of file scale_file_grads.F90.

863  integer, intent(in) :: fid, irecl
864  integer, intent(in) :: nx, ny, nz, k
865  logical, intent(in) :: yrev
866 
867  real(RP), intent(out) :: var(:)
868  integer, intent(out) :: ierr
869 
870  integer(1) :: buf(nx,ny)
871  integer :: i, j
872 
873  read(fid, rec=irecl, iostat=ierr) buf(:,:)
874  if ( ierr /= 0 ) return
875 
876  if ( yrev ) then
877  !$omp parallel do collapse(2)
878  do j = 1, ny
879  do i = 1, nx
880  var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,ny-j+1)
881  end do
882  end do
883  else
884  !$omp parallel do collapse(2)
885  do j = 1, ny
886  do i = 1, nx
887  var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,j)
888  end do
889  end do
890  end if
891 
892  return

References scale_prc::prc_abort().

Referenced by file_grads_close().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ upcase()

character(len=len(str)) function scale_file_grads::upcase ( character(len=*), intent(in)  str)

Definition at line 1093 of file scale_file_grads.F90.

1093  character(len=*), intent(in) :: str
1094  character(len=len(str)) :: upcase
1095  integer :: i
1096  do i = 1, len_trim(str)
1097  if ( str(i:i) >= 'a' .and. str(i:i) <= 'z' ) then
1098  upcase(i:i) = char( ichar(str(i:i)) - 32 )
1099  else
1100  upcase(i:i) = str(i:i)
1101  end if
1102  end do
1103  do i = len_trim(str)+1, len(str)
1104  upcase(i:i) = ' '
1105  end do

Referenced by file_grads_open(), and file_grads_varid().

Here is the caller graph for this function:
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_file_grads::file_grads_isoned
logical function, public file_grads_isoned(file_id, var_id)
Definition: scale_file_grads.F90:332
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11