MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
NCArrayReader.f90
Go to the documentation of this file.
1 !> @brief This module contains the NCArrayReaderModule
2 !!
3 !! This module defines the netcdf_array_load interface
4 !! which can read layered (UGRID) and non-layered (STRUCTURED)
5 !! netcdf arrays stored in modflow6 designated input variables.
6 !!
7 !<
9 
10  use kindmodule, only: dp, i4b, lgp
11  use constantsmodule, only: linelength
12  use simvariablesmodule, only: errmsg
18  use netcdfcommonmodule, only: nf_verify
19  use netcdf
20 
21  implicit none
22  private
23  public :: netcdf_array_load
24 
26  module procedure nc_array_load_int1d, nc_array_load_int2d, &
29  end interface netcdf_array_load
30 
31 contains
32 
33  !> @brief does the grid support per layer variables
34  !<
35  function is_layered(grid) result(layered)
36  character(len=*), intent(in) :: grid
37  logical(LGP) :: layered
38  select case (grid)
39  case ('LAYERED MESH')
40  layered = .true.
41  case ('STRUCTURED')
42  layered = .false.
43  case default
44  layered = .false.
45  end select
46  end function is_layered
47 
48  !> @brief Load NetCDF integer 1D array
49  !<
50  subroutine nc_array_load_int1d(int1d, mshape, idt, mf6_input, nc_vars, &
51  input_fname, iout, kper)
52  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: int1d
53  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape !< model shape
54  type(inputparamdefinitiontype), intent(in) :: idt !< input data type object describing this record
55  type(modflowinputtype), intent(in) :: mf6_input
56  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
57  character(len=*), intent(in) :: input_fname
58  integer(I4B), intent(in) :: iout
59  integer(I4B), optional, intent(in) :: kper
60  integer(I4B) :: varid
61  logical(LGP) :: layered
62 
63  layered = (idt%layered .and. is_layered(nc_vars%grid))
64 
65  if (layered) then
66  call load_integer1d_layered(int1d, mf6_input, mshape, idt, nc_vars, &
67  input_fname)
68  else
69  if (present(kper)) then
70  varid = nc_vars%varid(idt%mf6varname, period=kper)
71  else
72  varid = nc_vars%varid(idt%mf6varname)
73  end if
74  call load_integer1d_type(int1d, mf6_input, mshape, idt, nc_vars, &
75  varid, input_fname)
76  end if
77  end subroutine nc_array_load_int1d
78 
79  !> @brief Load NetCDF integer 2D array
80  !<
81  subroutine nc_array_load_int2d(int2d, mshape, idt, mf6_input, nc_vars, &
82  input_fname, iout)
83  integer(I4B), dimension(:, :), pointer, contiguous, intent(in) :: int2d
84  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape !< model shape
85  type(inputparamdefinitiontype), intent(in) :: idt !< input data type object describing this record
86  type(modflowinputtype), intent(in) :: mf6_input
87  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
88  character(len=*), intent(in) :: input_fname
89  integer(I4B), intent(in) :: iout
90  integer(I4B) :: varid
91  logical(LGP) :: layered
92 
93  layered = (idt%layered .and. is_layered(nc_vars%grid))
94 
95  if (layered) then
96  call load_integer2d_layered(int2d, mf6_input, mshape, idt, nc_vars, &
97  input_fname)
98  else
99  varid = nc_vars%varid(idt%mf6varname)
100  call load_integer2d_type(int2d, mf6_input, mshape, idt, nc_vars, &
101  varid, input_fname)
102  end if
103  end subroutine nc_array_load_int2d
104 
105  !> @brief Load NetCDF integer 3D array
106  !<
107  subroutine nc_array_load_int3d(int3d, mshape, idt, mf6_input, nc_vars, &
108  input_fname, iout)
109  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(in) :: int3d
110  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape !< model shape
111  type(inputparamdefinitiontype), intent(in) :: idt !< input data type object describing this record
112  type(modflowinputtype), intent(in) :: mf6_input
113  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
114  character(len=*), intent(in) :: input_fname
115  integer(I4B), intent(in) :: iout
116  integer(I4B) :: varid
117  logical(LGP) :: layered
118 
119  layered = (idt%layered .and. is_layered(nc_vars%grid))
120 
121  if (layered) then
122  call load_integer3d_layered(int3d, mf6_input, mshape, idt, nc_vars, &
123  input_fname)
124  else
125  varid = nc_vars%varid(idt%mf6varname)
126  call load_integer3d_type(int3d, mf6_input, mshape, idt, nc_vars, &
127  varid, input_fname)
128  end if
129  end subroutine nc_array_load_int3d
130 
131  !> @brief Load NetCDF double 1D array
132  !<
133  subroutine nc_array_load_dbl1d(dbl1d, mshape, idt, mf6_input, nc_vars, &
134  input_fname, iout, kper, iaux)
135  real(DP), dimension(:), pointer, contiguous, intent(in) :: dbl1d
136  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape !< model shape
137  type(inputparamdefinitiontype), intent(in) :: idt !< input data type object describing this record
138  type(modflowinputtype), intent(in) :: mf6_input
139  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
140  character(len=*), intent(in) :: input_fname
141  integer(I4B), intent(in) :: iout
142  integer(I4B), optional, intent(in) :: kper
143  integer(I4B), optional, intent(in) :: iaux
144  integer(I4B) :: varid
145  logical(LGP) :: layered
146 
147  if (present(kper)) then
148  layered = (kper > 0 .and. is_layered(nc_vars%grid))
149  else
150  layered = (idt%layered .and. is_layered(nc_vars%grid))
151  end if
152 
153  if (layered) then
154  if (present(kper)) then
155  call load_double1d_layered_spd(dbl1d, mf6_input, mshape, idt, nc_vars, &
156  kper, input_fname, iaux)
157  else
158  call load_double1d_layered(dbl1d, mf6_input, mshape, idt, nc_vars, &
159  input_fname)
160  end if
161  else
162  if (present(kper)) then
163  call load_double1d_spd(dbl1d, mf6_input, mshape, idt, nc_vars, &
164  kper, input_fname, iaux)
165  else
166  varid = nc_vars%varid(idt%mf6varname)
167  call load_double1d_type(dbl1d, mf6_input, mshape, idt, nc_vars, &
168  varid, input_fname)
169  end if
170  end if
171  end subroutine nc_array_load_dbl1d
172 
173  !> @brief Load NetCDF double 2D array
174  !<
175  subroutine nc_array_load_dbl2d(dbl2d, mshape, idt, mf6_input, nc_vars, &
176  input_fname, iout)
177  real(DP), dimension(:, :), pointer, contiguous, intent(in) :: dbl2d
178  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape !< model shape
179  type(inputparamdefinitiontype), intent(in) :: idt !< input data type object describing this record
180  type(modflowinputtype), intent(in) :: mf6_input
181  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
182  character(len=*), intent(in) :: input_fname
183  integer(I4B), intent(in) :: iout
184  integer(I4B) :: varid
185  logical(LGP) :: layered
186 
187  layered = (idt%layered .and. is_layered(nc_vars%grid))
188 
189  if (layered) then
190  call load_double2d_layered(dbl2d, mf6_input, mshape, idt, nc_vars, &
191  input_fname)
192  else
193  varid = nc_vars%varid(idt%mf6varname)
194  call load_double2d_type(dbl2d, mf6_input, mshape, idt, nc_vars, &
195  varid, input_fname)
196  end if
197  end subroutine nc_array_load_dbl2d
198 
199  !> @brief Load NetCDF double 3D array
200  !<
201  subroutine nc_array_load_dbl3d(dbl3d, mshape, idt, mf6_input, nc_vars, &
202  input_fname, iout)
203  real(DP), dimension(:, :, :), pointer, contiguous, intent(in) :: dbl3d
204  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape !< model shape
205  type(inputparamdefinitiontype), intent(in) :: idt !< input data type object describing this record
206  type(modflowinputtype), intent(in) :: mf6_input
207  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
208  character(len=*), intent(in) :: input_fname
209  integer(I4B), intent(in) :: iout
210  integer(I4B) :: varid
211  logical(LGP) :: layered
212 
213  layered = (idt%layered .and. is_layered(nc_vars%grid))
214 
215  if (layered) then
216  call load_double3d_layered(dbl3d, mf6_input, mshape, idt, nc_vars, &
217  input_fname)
218  else
219  varid = nc_vars%varid(idt%mf6varname)
220  call load_double3d_type(dbl3d, mf6_input, mshape, idt, nc_vars, &
221  varid, input_fname)
222  end if
223  end subroutine nc_array_load_dbl3d
224 
225  !> @brief load type 1d integer
226  !<
227  subroutine load_integer1d_type(int1d, mf6_input, mshape, idt, nc_vars, &
228  varid, input_fname)
229  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: int1d
230  type(modflowinputtype), intent(in) :: mf6_input
231  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
232  type(inputparamdefinitiontype), intent(in) :: idt
233  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
234  integer(I4B), intent(in) :: varid
235  character(len=*), intent(in) :: input_fname
236  integer(I4B), dimension(:), allocatable :: array_shape
237  integer(I4B), dimension(:, :, :), contiguous, pointer :: int3d_ptr
238  integer(I4B), dimension(:, :), contiguous, pointer :: int2d_ptr
239  integer(I4B) :: nvals
240 
241  ! initialize
242  nvals = 0
243 
244  if (idt%shape == 'NODES') then
245  ! set number of values
246  nvals = product(mshape)
247  if (size(mshape) == 3) then
248  int3d_ptr(1:mshape(3), 1:mshape(2), 1:mshape(1)) => int1d(1:nvals)
249  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int3d_ptr), &
250  nc_vars%nc_fname)
251  else if (size(mshape) == 2) then
252  int2d_ptr(1:mshape(2), 1:mshape(1)) => int1d(1:nvals)
253  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int2d_ptr), &
254  nc_vars%nc_fname)
255  else if (size(mshape) == 1) then
256  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d), nc_vars%nc_fname)
257  end if
258  else
259  ! interpret shape
260  call get_shape_from_string(idt%shape, array_shape, mf6_input%mempath)
261  ! set nvals
262  nvals = array_shape(1)
263  ! read and set data
264  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d), nc_vars%nc_fname)
265  end if
266  end subroutine load_integer1d_type
267 
268  !> @brief load type 1d integer layered
269  !<
270  subroutine load_integer1d_layered(int1d, mf6_input, mshape, idt, nc_vars, &
271  input_fname)
272  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: int1d
273  type(modflowinputtype), intent(in) :: mf6_input
274  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
275  type(inputparamdefinitiontype), intent(in) :: idt
276  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
277  character(len=*), intent(in) :: input_fname
278  integer(I4B), dimension(:), allocatable :: layer_shape
279  integer(I4B) :: nlay, varid
280  integer(I4B) :: k, ncpl
281  integer(I4B) :: index_start, index_stop
282  integer(I4B), dimension(:), contiguous, pointer :: int1d_ptr
283 
284  nullify (int1d_ptr)
285 
286  call get_layered_shape(mshape, nlay, layer_shape)
287 
288  ncpl = product(layer_shape)
289  index_start = 1
290  do k = 1, nlay
291  varid = nc_vars%varid(idt%mf6varname, layer=k)
292  index_stop = index_start + ncpl - 1
293  int1d_ptr(1:ncpl) => int1d(index_start:index_stop)
294  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
295  nc_vars%nc_fname)
296  index_start = index_stop + 1
297  end do
298  end subroutine load_integer1d_layered
299 
300  !> @brief load type 2d integer
301  !<
302  subroutine load_integer2d_type(int2d, mf6_input, mshape, idt, nc_vars, varid, &
303  input_fname)
304  integer(I4B), dimension(:, :), contiguous, pointer, intent(in) :: int2d
305  type(modflowinputtype), intent(in) :: mf6_input
306  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
307  type(inputparamdefinitiontype), intent(in) :: idt
308  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
309  integer(I4B), intent(in) :: varid
310  character(len=*), intent(in) :: input_fname
311  integer(I4B), dimension(:), contiguous, pointer :: int1d_ptr
312  integer(I4B), dimension(:), allocatable :: array_shape
313  integer(I4B) :: ncpl, nlay
314 
315  nullify (int1d_ptr)
316 
317  if (nc_vars%grid == 'STRUCTURED') then
318  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int2d), nc_vars%nc_fname)
319  else if (nc_vars%grid == 'LAYERED MESH') then
320  call get_layered_shape(mshape, nlay, array_shape)
321  ncpl = product(array_shape)
322  int1d_ptr(1:ncpl) => int2d(:, :)
323  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
324  nc_vars%nc_fname)
325  end if
326  end subroutine load_integer2d_type
327 
328  !> @brief load type 2d integer layered
329  !<
330  subroutine load_integer2d_layered(int2d, mf6_input, mshape, idt, nc_vars, &
331  input_fname)
332  integer(I4B), dimension(:, :), contiguous, pointer, intent(in) :: int2d
333  type(modflowinputtype), intent(in) :: mf6_input
334  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
335  type(inputparamdefinitiontype), intent(in) :: idt
336  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
337  character(len=*), intent(in) :: input_fname
338  integer(I4B), dimension(:), allocatable :: layer_shape
339  integer(I4B) :: k
340  integer(I4B) :: ncpl, nlay, varid
341  integer(I4B), dimension(:), contiguous, pointer :: int1d_ptr
342 
343  nullify (int1d_ptr)
344 
345  if (size(mshape) == 3) then
346  write (errmsg, '(a,a,a)') &
347  'Layered netcdf read not supported for DIS int2d type ('// &
348  trim(idt%tagname)//').'
349  call store_error(errmsg)
350  call store_error_filename(input_fname)
351  else if (size(mshape) == 2) then
352  call get_layered_shape(mshape, nlay, layer_shape)
353  ncpl = layer_shape(1)
354  do k = 1, nlay
355  varid = nc_vars%varid(idt%mf6varname, layer=k)
356  int1d_ptr(1:ncpl) => int2d(1:ncpl, k)
357  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
358  nc_vars%nc_fname)
359  end do
360  end if
361  end subroutine load_integer2d_layered
362 
363  !> @brief load type 3d integer
364  !<
365  subroutine load_integer3d_type(int3d, mf6_input, mshape, idt, nc_vars, varid, &
366  input_fname)
367  integer(I4B), dimension(:, :, :), contiguous, pointer, intent(in) :: int3d
368  type(modflowinputtype), intent(in) :: mf6_input
369  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
370  type(inputparamdefinitiontype), intent(in) :: idt
371  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
372  integer(I4B), intent(in) :: varid
373  character(len=*), intent(in) :: input_fname
374  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int3d), nc_vars%nc_fname)
375  end subroutine load_integer3d_type
376 
377  !> @brief load type 3d integer layered
378  !<
379  subroutine load_integer3d_layered(int3d, mf6_input, mshape, idt, nc_vars, &
380  input_fname)
381  integer(I4B), dimension(:, :, :), contiguous, pointer, intent(in) :: int3d
382  type(modflowinputtype), intent(in) :: mf6_input
383  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
384  type(inputparamdefinitiontype), intent(in) :: idt
385  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
386  character(len=*), intent(in) :: input_fname
387  integer(I4B), dimension(:), allocatable :: layer_shape
388  integer(I4B) :: k !, i, j
389  integer(I4B) :: ncpl, nlay, varid
390  integer(I4B) :: index_start, index_stop
391  integer(I4B), dimension(:), contiguous, pointer :: int1d_ptr
392 
393  nullify (int1d_ptr)
394  index_start = 1
395  call get_layered_shape(mshape, nlay, layer_shape)
396  ncpl = product(layer_shape)
397 
398  do k = 1, nlay
399  varid = nc_vars%varid(idt%mf6varname, layer=k)
400  index_stop = index_start + ncpl - 1
401  int1d_ptr(1:ncpl) => int3d(:, :, k:k)
402  call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
403  nc_vars%nc_fname)
404  index_start = index_stop + 1
405  end do
406  end subroutine load_integer3d_layered
407 
408  !> @brief load type 1d double
409  !<
410  subroutine load_double1d_type(dbl1d, mf6_input, mshape, idt, nc_vars, &
411  varid, input_fname)
412  real(DP), dimension(:), contiguous, pointer, intent(in) :: dbl1d
413  type(modflowinputtype), intent(in) :: mf6_input
414  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
415  type(inputparamdefinitiontype), intent(in) :: idt
416  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
417  integer(I4B), intent(in) :: varid
418  character(len=*), intent(in) :: input_fname
419  integer(I4B), dimension(:), allocatable :: array_shape
420  real(DP), dimension(:, :, :), contiguous, pointer :: dbl3d_ptr
421  real(DP), dimension(:, :), contiguous, pointer :: dbl2d_ptr
422  integer(I4B) :: nvals
423 
424  ! initialize
425  nvals = 0
426 
427  if (idt%shape == 'NODES') then
428  ! set number of values
429  nvals = product(mshape)
430  if (size(mshape) == 3) then
431  dbl3d_ptr(1:mshape(3), 1:mshape(2), 1:mshape(1)) => dbl1d(1:nvals)
432  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl3d_ptr), &
433  nc_vars%nc_fname)
434  else if (size(mshape) == 2) then
435  dbl2d_ptr(1:mshape(2), 1:mshape(1)) => dbl1d(1:nvals)
436  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl2d_ptr), &
437  nc_vars%nc_fname)
438  else if (size(mshape) == 1) then
439  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d), nc_vars%nc_fname)
440  end if
441  else
442  ! interpret shape
443  call get_shape_from_string(idt%shape, array_shape, mf6_input%mempath)
444  ! set nvals
445  nvals = array_shape(1)
446  ! read and set data
447  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d), nc_vars%nc_fname)
448  end if
449  end subroutine load_double1d_type
450 
451  !> @brief load type 1d double
452  !<
453  subroutine load_double1d_spd(dbl1d, mf6_input, mshape, idt, nc_vars, &
454  iper, input_fname, iaux)
455  use constantsmodule, only: dnodata
456  real(DP), dimension(:), contiguous, pointer, intent(in) :: dbl1d
457  type(modflowinputtype), intent(in) :: mf6_input
458  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
459  type(inputparamdefinitiontype), intent(in) :: idt
460  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
461  integer(I4B), intent(in) :: iper
462  character(len=*), intent(in) :: input_fname
463  integer(I4B), optional, intent(in) :: iaux
464  real(DP), dimension(:, :, :), contiguous, pointer :: dbl3d
465  integer(I4B) :: nvals, varid
466  integer(I4B) :: n, i, j, k
467 
468  ! initialize
469  nvals = 0
470 
471  ! set varid
472  if (present(iaux)) then
473  varid = nc_vars%varid(idt%mf6varname, period=iper, iaux=iaux)
474  else
475  varid = nc_vars%varid(idt%mf6varname, period=iper)
476  end if
477 
478  if (idt%shape == 'NODES') then
479  ! TODO future support
480  write (errmsg, '(a)') &
481  'IDM NetCDF load_double1d_spd NODES var shape not supported => '// &
482  trim(idt%tagname)
483  call store_error(errmsg)
484  call store_error_filename(input_fname)
485  else if (idt%shape == 'NCPL' .or. idt%shape == 'NAUX NCPL') then
486 
487  if (size(mshape) == 3) then
488  allocate (dbl3d(mshape(3), mshape(2), mshape(1)))
489  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl3d), &
490  nc_vars%nc_fname)
491  n = 0
492  do k = 1, size(dbl3d, dim=3)
493  do i = 1, size(dbl3d, dim=2)
494  do j = 1, size(dbl3d, dim=1)
495  if (n < size(dbl1d)) then
496  n = n + 1
497  else
498  n = 1
499  end if
500  if (dbl3d(j, i, k) /= dnodata) then
501  dbl1d(n) = dbl3d(j, i, k)
502  end if
503  end do
504  end do
505  end do
506 
507  else if (size(mshape) == 2) then
508  ! TODO
509  write (errmsg, '(a)') &
510  'IDM NetCDF load_double1d_spd DISV model not supported => '// &
511  trim(idt%tagname)
512  call store_error(errmsg)
513  call store_error_filename(input_fname)
514  end if
515  end if
516  end subroutine load_double1d_spd
517 
518  !> @brief load type 1d double layered
519  !<
520  subroutine load_double1d_layered(dbl1d, mf6_input, mshape, idt, nc_vars, &
521  input_fname)
522  real(DP), dimension(:), contiguous, pointer, intent(in) :: dbl1d
523  type(modflowinputtype), intent(in) :: mf6_input
524  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
525  type(inputparamdefinitiontype), intent(in) :: idt
526  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
527  character(len=*), intent(in) :: input_fname
528  integer(I4B), dimension(:), allocatable :: layer_shape
529  integer(I4B) :: nlay, varid
530  integer(I4B) :: k, ncpl
531  integer(I4B) :: index_start, index_stop
532  real(DP), dimension(:), contiguous, pointer :: dbl1d_ptr
533 
534  nullify (dbl1d_ptr)
535  index_start = 1
536  call get_layered_shape(mshape, nlay, layer_shape)
537  ncpl = product(layer_shape)
538 
539  do k = 1, nlay
540  varid = nc_vars%varid(idt%mf6varname, layer=k)
541  index_stop = index_start + ncpl - 1
542  dbl1d_ptr(1:ncpl) => dbl1d(index_start:index_stop)
543  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
544  nc_vars%nc_fname)
545  index_start = index_stop + 1
546  end do
547  end subroutine load_double1d_layered
548 
549  !> @brief load type 1d double layered
550  !<
551  subroutine load_double1d_layered_spd(dbl1d, mf6_input, mshape, idt, nc_vars, &
552  iper, input_fname, iaux)
553  use constantsmodule, only: dnodata
554  real(DP), dimension(:), contiguous, pointer, intent(in) :: dbl1d
555  type(modflowinputtype), intent(in) :: mf6_input
556  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
557  type(inputparamdefinitiontype), intent(in) :: idt
558  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
559  integer(I4B), intent(in) :: iper
560  character(len=*), intent(in) :: input_fname
561  integer(I4B), optional, intent(in) :: iaux
562  integer(I4B), dimension(:), allocatable :: layer_shape
563  integer(I4B) :: nlay, varid
564  integer(I4B) :: k, n, ncpl
565  integer(I4B) :: index_start, index_stop
566  real(DP), dimension(:), contiguous, pointer :: dbl1d_ptr
567 
568  call get_layered_shape(mshape, nlay, layer_shape)
569  ncpl = product(layer_shape)
570  allocate (dbl1d_ptr(ncpl))
571 
572  do k = 1, nlay
573  index_start = 1
574  index_stop = index_start + ncpl - 1
575  if (present(iaux)) then
576  varid = nc_vars%varid(idt%mf6varname, layer=k, period=iper, iaux=iaux)
577  else
578  varid = nc_vars%varid(idt%mf6varname, layer=k, period=iper)
579  end if
580  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
581  nc_vars%nc_fname)
582  do n = 1, ncpl
583  if (dbl1d_ptr(n) /= dnodata) then
584  dbl1d(n) = dbl1d_ptr(n)
585  end if
586  end do
587  end do
588 
589  ! cleanup
590  deallocate (dbl1d_ptr)
591  end subroutine load_double1d_layered_spd
592 
593  !> @brief load type 2d double
594  !<
595  subroutine load_double2d_type(dbl2d, mf6_input, mshape, idt, nc_vars, varid, &
596  input_fname)
597  real(DP), dimension(:, :), contiguous, pointer, intent(in) :: dbl2d
598  type(modflowinputtype), intent(in) :: mf6_input
599  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
600  type(inputparamdefinitiontype), intent(in) :: idt
601  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
602  integer(I4B), intent(in) :: varid
603  character(len=*), intent(in) :: input_fname
604  real(DP), dimension(:), contiguous, pointer :: dbl1d_ptr
605  integer(I4B), dimension(:), allocatable :: array_shape
606  integer(I4B) :: ncpl, nlay
607 
608  nullify (dbl1d_ptr)
609 
610  if (nc_vars%grid == 'STRUCTURED') then
611  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl2d), nc_vars%nc_fname)
612  else if (nc_vars%grid == 'LAYERED MESH') then
613  call get_layered_shape(mshape, nlay, array_shape)
614  ncpl = product(array_shape)
615  dbl1d_ptr(1:ncpl) => dbl2d(:, :)
616  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
617  nc_vars%nc_fname)
618  end if
619  end subroutine load_double2d_type
620 
621  !> @brief load type 2d double layered
622  !<
623  subroutine load_double2d_layered(dbl2d, mf6_input, mshape, idt, nc_vars, &
624  input_fname)
625  real(DP), dimension(:, :), contiguous, pointer, intent(in) :: dbl2d
626  type(modflowinputtype), intent(in) :: mf6_input
627  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
628  type(inputparamdefinitiontype), intent(in) :: idt
629  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
630  character(len=*), intent(in) :: input_fname
631  integer(I4B), dimension(:), allocatable :: layer_shape
632  integer(I4B) :: k
633  integer(I4B) :: ncpl, nlay, varid
634  real(DP), dimension(:), contiguous, pointer :: dbl1d_ptr
635 
636  nullify (dbl1d_ptr)
637 
638  if (size(mshape) == 3) then
639  write (errmsg, '(a,a,a)') &
640  'Layered netcdf read not supported for DIS dbl2d type ('// &
641  trim(idt%tagname)//').'
642  call store_error(errmsg)
643  call store_error_filename(input_fname)
644  else if (size(mshape) == 2) then
645  call get_layered_shape(mshape, nlay, layer_shape)
646  ncpl = layer_shape(1)
647  do k = 1, nlay
648  varid = nc_vars%varid(idt%mf6varname, layer=k)
649  dbl1d_ptr(1:ncpl) => dbl2d(1:ncpl, k)
650  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
651  nc_vars%nc_fname)
652  end do
653  end if
654  end subroutine load_double2d_layered
655 
656  !> @brief load type 3d double
657  !<
658  subroutine load_double3d_type(dbl3d, mf6_input, mshape, idt, nc_vars, varid, &
659  input_fname)
660  real(DP), dimension(:, :, :), contiguous, pointer, intent(in) :: dbl3d
661  type(modflowinputtype), intent(in) :: mf6_input
662  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
663  type(inputparamdefinitiontype), intent(in) :: idt
664  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
665  integer(I4B), intent(in) :: varid
666  character(len=*), intent(in) :: input_fname
667  !
668  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl3d), nc_vars%nc_fname)
669  end subroutine load_double3d_type
670 
671  !> @brief load type 3d double layered
672  !<
673  subroutine load_double3d_layered(dbl3d, mf6_input, mshape, idt, nc_vars, &
674  input_fname)
675  real(DP), dimension(:, :, :), contiguous, pointer, intent(in) :: dbl3d
676  type(modflowinputtype), intent(in) :: mf6_input
677  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: mshape
678  type(inputparamdefinitiontype), intent(in) :: idt
679  type(ncpackagevarstype), pointer, intent(in) :: nc_vars
680  character(len=*), intent(in) :: input_fname
681  integer(I4B), dimension(:), allocatable :: layer_shape
682  integer(I4B) :: k !, i, j
683  integer(I4B) :: ncpl, nlay, varid
684  integer(I4B) :: index_start, index_stop
685  real(DP), dimension(:), contiguous, pointer :: dbl1d_ptr
686 
687  nullify (dbl1d_ptr)
688 
689  call get_layered_shape(mshape, nlay, layer_shape)
690 
691  ncpl = product(layer_shape)
692  index_start = 1
693  do k = 1, nlay
694  varid = nc_vars%varid(idt%mf6varname, layer=k)
695  index_stop = index_start + ncpl - 1
696  dbl1d_ptr(1:ncpl) => dbl3d(:, :, k:k)
697  call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
698  nc_vars%nc_fname)
699  index_start = index_stop + 1
700  end do
701  end subroutine load_double3d_layered
702 
703 end module ncarrayreadermodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
This module contains the InputDefinitionModule.
This module defines variable data types.
Definition: kind.f90:8
This module contains the ModflowInputModule.
Definition: ModflowInput.f90:9
This module contains the NCArrayReaderModule.
subroutine load_double1d_layered(dbl1d, mf6_input, mshape, idt, nc_vars, input_fname)
load type 1d double layered
subroutine nc_array_load_int2d(int2d, mshape, idt, mf6_input, nc_vars, input_fname, iout)
Load NetCDF integer 2D array.
subroutine load_double3d_layered(dbl3d, mf6_input, mshape, idt, nc_vars, input_fname)
load type 3d double layered
subroutine nc_array_load_dbl2d(dbl2d, mshape, idt, mf6_input, nc_vars, input_fname, iout)
Load NetCDF double 2D array.
subroutine load_integer3d_layered(int3d, mf6_input, mshape, idt, nc_vars, input_fname)
load type 3d integer layered
subroutine load_integer1d_type(int1d, mf6_input, mshape, idt, nc_vars, varid, input_fname)
load type 1d integer
subroutine load_double3d_type(dbl3d, mf6_input, mshape, idt, nc_vars, varid, input_fname)
load type 3d double
logical(lgp) function is_layered(grid)
does the grid support per layer variables
subroutine nc_array_load_int3d(int3d, mshape, idt, mf6_input, nc_vars, input_fname, iout)
Load NetCDF integer 3D array.
subroutine load_double2d_type(dbl2d, mf6_input, mshape, idt, nc_vars, varid, input_fname)
load type 2d double
subroutine load_integer3d_type(int3d, mf6_input, mshape, idt, nc_vars, varid, input_fname)
load type 3d integer
subroutine load_integer2d_layered(int2d, mf6_input, mshape, idt, nc_vars, input_fname)
load type 2d integer layered
subroutine load_integer2d_type(int2d, mf6_input, mshape, idt, nc_vars, varid, input_fname)
load type 2d integer
subroutine nc_array_load_int1d(int1d, mshape, idt, mf6_input, nc_vars, input_fname, iout, kper)
Load NetCDF integer 1D array.
subroutine load_double1d_layered_spd(dbl1d, mf6_input, mshape, idt, nc_vars, iper, input_fname, iaux)
load type 1d double layered
subroutine load_integer1d_layered(int1d, mf6_input, mshape, idt, nc_vars, input_fname)
load type 1d integer layered
subroutine load_double1d_type(dbl1d, mf6_input, mshape, idt, nc_vars, varid, input_fname)
load type 1d double
subroutine load_double2d_layered(dbl2d, mf6_input, mshape, idt, nc_vars, input_fname)
load type 2d double layered
subroutine nc_array_load_dbl3d(dbl3d, mshape, idt, mf6_input, nc_vars, input_fname, iout)
Load NetCDF double 3D array.
subroutine nc_array_load_dbl1d(dbl1d, mshape, idt, mf6_input, nc_vars, input_fname, iout, kper, iaux)
Load NetCDF double 1D array.
subroutine load_double1d_spd(dbl1d, mf6_input, mshape, idt, nc_vars, iper, input_fname, iaux)
load type 1d double
This module contains the NCFileVarsModule.
Definition: NCFileVars.f90:7
This module contains the NetCDFCommonModule.
Definition: NetCDFCommon.f90:6
subroutine, public nf_verify(res, nc_fname)
error check a netcdf-fortran interface call
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
subroutine, public get_layered_shape(mshape, nlay, layer_shape)
subroutine, public get_shape_from_string(shape_string, array_shape, memoryPath)
derived type for storing input definition for a file
Type describing input variables for a package in NetCDF file.
Definition: NCFileVars.f90:22