MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
gridfilereadermodule Module Reference

Data Types

type  gridfilereadertype
 

Functions/Subroutines

subroutine initialize (this, iu)
 @Brief Initialize the grid file reader. More...
 
subroutine finalize (this)
 Finalize the grid file reader. More...
 
subroutine read_header (this)
 Read the file's self-describing header. Internal use only. More...
 
subroutine read_header_meta (this)
 Read self-describing metadata (first four lines). Internal use only. More...
 
subroutine read_header_body (this)
 Read the header body section (text following first. More...
 
integer(i4b) function read_int (this, key)
 Read an integer scalar from a grid file. More...
 
real(dp) function read_dbl (this, key)
 Read a double precision scalar from a grid file. More...
 
integer(i4b) function, dimension(:), allocatable read_int_1d (this, key)
 Read a 1D integer array from a grid file. More...
 
real(dp) function, dimension(:), allocatable read_dbl_1d (this, key)
 Read a 1D double array from a grid file. More...
 
character(len=:) function, allocatable read_charstr (this, key)
 Read a character string from a grid file. More...
 
integer(i4b) function, dimension(:), allocatable read_grid_shape (this)
 Read the grid shape from a grid file. More...
 

Function/Subroutine Documentation

◆ finalize()

subroutine gridfilereadermodule::finalize ( class(gridfilereadertype this)

Definition at line 62 of file GridFileReader.f90.

63  class(GridFileReaderType) :: this
64 
65  close (this%inunit)
66  call hash_table_da(this%dim)
67  call hash_table_da(this%pos)
68  call hash_table_da(this%typ)
69  call hash_table_da(this%shp_idx)
70  deallocate (this%shp)
71 
Here is the call graph for this function:

◆ initialize()

subroutine gridfilereadermodule::initialize ( class(gridfilereadertype this,
integer(i4b), intent(in)  iu 
)

Definition at line 47 of file GridFileReader.f90.

48  class(GridFileReaderType) :: this
49  integer(I4B), intent(in) :: iu
50 
51  this%inunit = iu
52  call hash_table_cr(this%dim)
53  call hash_table_cr(this%pos)
54  call hash_table_cr(this%typ)
55  call hash_table_cr(this%shp_idx)
56  allocate (this%shp(0))
57  call this%read_header()
58 
Here is the call graph for this function:

◆ read_charstr()

character(len=:) function, allocatable gridfilereadermodule::read_charstr ( class(gridfilereadertype), intent(inout)  this,
character(len=*), intent(in)  key 
)

Definition at line 323 of file GridFileReader.f90.

324  class(GridFileReaderType), intent(inout) :: this
325  character(len=*), intent(in) :: key
326  character(len=:), allocatable :: charstr
327  ! local
328  integer(I4B) :: idx, ndim, nvals, pos, typ
329  character(len=:), allocatable :: msg
330 
331  msg = 'Variable '//trim(key)//' is not a character array'
332  ndim = this%dim%get(key)
333  if (ndim /= 1) then
334  write (errmsg, '(a)') msg
335  call store_error(errmsg, terminate=.true.)
336  end if
337  typ = this%typ%get(key)
338  if (typ /= 3) then
339  write (errmsg, '(a)') msg
340  call store_error(errmsg, terminate=.true.)
341  end if
342  idx = this%shp_idx%get(key)
343  pos = this%pos%get(key)
344  nvals = this%shp(idx)
345  allocate (character(nvals) :: charstr)
346  read (this%inunit, pos=pos) charstr
347  rewind(this%inunit)
Here is the call graph for this function:

◆ read_dbl()

real(dp) function gridfilereadermodule::read_dbl ( class(gridfilereadertype), intent(inout)  this,
character(len=*), intent(in)  key 
)

Definition at line 239 of file GridFileReader.f90.

240  class(GridFileReaderType), intent(inout) :: this
241  character(len=*), intent(in) :: key
242  real(DP) :: v
243  ! local
244  integer(I4B) :: ndim, pos, typ
245  character(len=:), allocatable :: msg
246 
247  msg = 'Variable '//trim(key)//' is not a double precision scalar'
248  ndim = this%dim%get(key)
249  if (ndim /= 0) then
250  write (errmsg, '(a)') msg
251  call store_error(errmsg, terminate=.true.)
252  end if
253  typ = this%typ%get(key)
254  if (typ /= 2) then
255  write (errmsg, '(a)') msg
256  call store_error(errmsg, terminate=.true.)
257  end if
258  pos = this%pos%get(key)
259  read (this%inunit, pos=pos) v
260  rewind(this%inunit)
261 
Here is the call graph for this function:

◆ read_dbl_1d()

real(dp) function, dimension(:), allocatable gridfilereadermodule::read_dbl_1d ( class(gridfilereadertype), intent(inout)  this,
character(len=*), intent(in)  key 
)

Definition at line 294 of file GridFileReader.f90.

295  class(GridFileReaderType), intent(inout) :: this
296  character(len=*), intent(in) :: key
297  real(DP), allocatable :: v(:)
298  ! local
299  integer(I4B) :: idx, ndim, nvals, pos, typ
300  character(len=:), allocatable :: msg
301 
302  msg = 'Variable '//trim(key)//' is not a 1D double array'
303  ndim = this%dim%get(key)
304  if (ndim /= 1) then
305  write (errmsg, '(a)') msg
306  call store_error(errmsg, terminate=.true.)
307  end if
308  typ = this%typ%get(key)
309  if (typ /= 2) then
310  write (errmsg, '(a)') msg
311  call store_error(errmsg, terminate=.true.)
312  end if
313  idx = this%shp_idx%get(key)
314  pos = this%pos%get(key)
315  nvals = this%shp(idx)
316  allocate (v(nvals))
317  read (this%inunit, pos=pos) v
318  rewind(this%inunit)
319 
Here is the call graph for this function:

◆ read_grid_shape()

integer(i4b) function, dimension(:), allocatable gridfilereadermodule::read_grid_shape ( class(gridfilereadertype this)

Definition at line 351 of file GridFileReader.f90.

352  ! dummy
353  class(GridFileReaderType) :: this
354  integer(I4B), allocatable :: v(:)
355 
356  select case (this%grid_type)
357  case ("DIS")
358  allocate (v(3))
359  v(1) = this%read_int("NLAY")
360  v(2) = this%read_int("NROW")
361  v(3) = this%read_int("NCOL")
362  case ("DISV")
363  allocate (v(2))
364  v(1) = this%read_int("NLAY")
365  v(2) = this%read_int("NCPL")
366  case ("DISU")
367  allocate (v(1))
368  v(1) = this%read_int("NODES")
369  case ("DIS2D")
370  allocate (v(2))
371  v(1) = this%read_int("NROW")
372  v(2) = this%read_int("NCOL")
373  case ("DISV2D")
374  allocate (v(1))
375  v(1) = this%read_int("NODES")
376  case ("DISV1D")
377  allocate (v(1))
378  v(1) = this%read_int("NCELLS")
379  end select
380 

◆ read_header()

subroutine gridfilereadermodule::read_header ( class(gridfilereadertype this)

Definition at line 75 of file GridFileReader.f90.

76  class(GridFileReaderType) :: this
77  call this%read_header_meta()
78  call this%read_header_body()

◆ read_header_body()

subroutine gridfilereadermodule::read_header_body ( class(gridfilereadertype this)

Definition at line 128 of file GridFileReader.f90.

129  ! dummy
130  class(GridFileReaderType) :: this
131  ! local
132  character(len=:), allocatable :: body
133  character(len=:), allocatable :: line
134  character(len=10) :: key, dtype
135  real(DP) :: rval
136  integer(I4B) :: i, lloc, istart, istop, ival, pos
137  integer(I4B) :: nvars, ndim, dim, ishp
138  integer(I4B), allocatable :: shp(:)
139 
140  allocate (this%keys(this%ntxt))
141  allocate (character(len=this%lentxt*this%ntxt) :: body)
142  allocate (character(len=this%lentxt) :: line)
143 
144  nvars = 0
145  read (this%inunit) body
146  inquire (this%inunit, pos=pos)
147  do i = 1, this%lentxt * this%ntxt, this%lentxt
148  line = body(i:i + this%lentxt - 1)
149  lloc = 1
150 
151  ! key
152  lloc = 1
153  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
154  key = line(istart:istop)
155  nvars = nvars + 1
156  this%keys(nvars) = key
157 
158  ! type
159  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
160  dtype = line(istart:istop)
161  if (dtype == "INTEGER") then
162  call this%typ%add(key, 1)
163  else if (dtype == "DOUBLE") then
164  call this%typ%add(key, 2)
165  else if (dtype == "CHARACTER") then
166  call this%typ%add(key, 3)
167  end if
168 
169  ! dims
170  call urword(line, lloc, istart, istop, 0, ival, rval, 0, 0)
171  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
172  ndim = ival
173  call this%dim%add(key, ndim)
174 
175  ! shape
176  if (allocated(shp)) deallocate (shp)
177  allocate (shp(ndim))
178  if (ndim > 0) then
179  do dim = 1, ndim
180  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
181  shp(dim) = ival
182  end do
183  ishp = size(this%shp)
184  call expandarray(this%shp, increment=ndim)
185  this%shp(ishp + 1:ishp + ndim) = shp
186  call this%shp_idx%add(key, ishp + 1)
187  end if
188 
189  ! position
190  call this%pos%add(key, pos)
191  if (ndim == 0) then
192  if (dtype == "INTEGER") then
193  pos = pos + 4
194  else if (dtype == "DOUBLE") then
195  pos = pos + 8
196  end if
197  else
198  if (dtype == "INTEGER") then
199  pos = pos + (product(shp) * 4)
200  else if (dtype == "DOUBLE") then
201  pos = pos + (product(shp) * 8)
202  else if (dtype == "CHARACTER") then
203  pos = pos + (product(shp) * 8)
204  end if
205  end if
206  end do
207 
208  rewind(this%inunit)
209 
Here is the call graph for this function:

◆ read_header_meta()

subroutine gridfilereadermodule::read_header_meta ( class(gridfilereadertype this)

Definition at line 82 of file GridFileReader.f90.

83  ! dummy
84  class(GridFileReaderType) :: this
85  ! local
86  character(len=50) :: line
87  integer(I4B) :: lloc, istart, istop
88  integer(I4B) :: ival
89  real(DP) :: rval
90 
91  ! grid type
92  read (this%inunit) line
93  lloc = 1
94  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
95  if (line(istart:istop) /= 'GRID') then
96  call store_error('Binary grid file must begin with "GRID". '//&
97  &'Found: '//line(istart:istop))
98  call store_error_unit(this%inunit)
99  end if
100  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
101  this%grid_type = line(istart:istop)
102 
103  ! version
104  read (this%inunit) line
105  lloc = 1
106  call urword(line, lloc, istart, istop, 0, ival, rval, 0, 0)
107  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
108  this%version = ival
109 
110  ! ntxt
111  read (this%inunit) line
112  lloc = 1
113  call urword(line, lloc, istart, istop, 0, ival, rval, 0, 0)
114  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
115  this%ntxt = ival
116 
117  ! lentxt
118  read (this%inunit) line
119  lloc = 1
120  call urword(line, lloc, istart, istop, 0, ival, rval, 0, 0)
121  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
122  this%lentxt = ival
123 
Here is the call graph for this function:

◆ read_int()

integer(i4b) function gridfilereadermodule::read_int ( class(gridfilereadertype), intent(inout)  this,
character(len=*), intent(in)  key 
)

Definition at line 213 of file GridFileReader.f90.

214  class(GridFileReaderType), intent(inout) :: this
215  character(len=*), intent(in) :: key
216  integer(I4B) :: v
217  ! local
218  integer(I4B) :: ndim, pos, typ
219  character(len=:), allocatable :: msg
220 
221  msg = 'Variable '//trim(key)//' is not an integer scalar'
222  ndim = this%dim%get(key)
223  if (ndim /= 0) then
224  write (errmsg, '(a)') msg
225  call store_error(errmsg, terminate=.true.)
226  end if
227  typ = this%typ%get(key)
228  if (typ /= 1) then
229  write (errmsg, '(a)') msg
230  call store_error(errmsg, terminate=.true.)
231  end if
232  pos = this%pos%get(key)
233  read (this%inunit, pos=pos) v
234  rewind(this%inunit)
235 
Here is the call graph for this function:

◆ read_int_1d()

integer(i4b) function, dimension(:), allocatable gridfilereadermodule::read_int_1d ( class(gridfilereadertype), intent(inout)  this,
character(len=*), intent(in)  key 
)

Definition at line 265 of file GridFileReader.f90.

266  class(GridFileReaderType), intent(inout) :: this
267  character(len=*), intent(in) :: key
268  integer(I4B), allocatable :: v(:)
269  ! local
270  integer(I4B) :: idx, ndim, nvals, pos, typ
271  character(len=:), allocatable :: msg
272 
273  msg = 'Variable '//trim(key)//' is not a 1D integer array'
274  ndim = this%dim%get(key)
275  if (ndim /= 1) then
276  write (errmsg, '(a)') msg
277  call store_error(errmsg, terminate=.true.)
278  end if
279  typ = this%typ%get(key)
280  if (typ /= 1) then
281  write (errmsg, '(a)') msg
282  call store_error(errmsg, terminate=.true.)
283  end if
284  idx = this%shp_idx%get(key)
285  pos = this%pos%get(key)
286  nvals = this%shp(idx)
287  allocate (v(nvals))
288  read (this%inunit, pos=pos) v
289  rewind(this%inunit)
290 
Here is the call graph for this function: