36 character(len=*),
intent(in) :: grid
37 logical(LGP) :: layered
51 input_fname, iout, kper)
52 integer(I4B),
dimension(:),
pointer,
contiguous,
intent(in) :: int1d
53 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
57 character(len=*),
intent(in) :: input_fname
58 integer(I4B),
intent(in) :: iout
59 integer(I4B),
optional,
intent(in) :: kper
60 integer(I4B) :: varid, iper
61 logical(LGP) :: layered
64 layered = (idt%layered .and.
is_layered(nc_vars%grid))
66 if (
present(kper))
then
83 varid = nc_vars%varid(idt%tagname)
94 integer(I4B),
dimension(:, :),
pointer,
contiguous,
intent(in) :: int2d
95 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
99 character(len=*),
intent(in) :: input_fname
100 integer(I4B),
intent(in) :: iout
101 integer(I4B) :: varid
102 logical(LGP) :: layered
104 layered = (idt%layered .and.
is_layered(nc_vars%grid))
110 varid = nc_vars%varid(idt%tagname)
120 integer(I4B),
dimension(:, :, :),
pointer,
contiguous,
intent(in) :: int3d
121 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
125 character(len=*),
intent(in) :: input_fname
126 integer(I4B),
intent(in) :: iout
127 integer(I4B) :: varid
128 logical(LGP) :: layered
130 layered = (idt%layered .and.
is_layered(nc_vars%grid))
136 varid = nc_vars%varid(idt%tagname)
145 input_fname, iout, kper, iaux)
146 real(DP),
dimension(:),
pointer,
contiguous,
intent(in) :: dbl1d
147 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
151 character(len=*),
intent(in) :: input_fname
152 integer(I4B),
intent(in) :: iout
153 integer(I4B),
optional,
intent(in) :: kper
154 integer(I4B),
optional,
intent(in) :: iaux
155 integer(I4B) :: varid, iper
156 logical(LGP) :: layered
159 layered = (idt%layered .and.
is_layered(nc_vars%grid))
161 if (
present(kper))
then
168 iper, input_fname, iaux)
176 iper, input_fname, iaux)
178 varid = nc_vars%varid(idt%tagname)
189 real(DP),
dimension(:, :),
pointer,
contiguous,
intent(in) :: dbl2d
190 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
194 character(len=*),
intent(in) :: input_fname
195 integer(I4B),
intent(in) :: iout
196 integer(I4B) :: varid
197 logical(LGP) :: layered
199 layered = (idt%layered .and.
is_layered(nc_vars%grid))
205 varid = nc_vars%varid(idt%tagname)
215 real(DP),
dimension(:, :, :),
pointer,
contiguous,
intent(in) :: dbl3d
216 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
220 character(len=*),
intent(in) :: input_fname
221 integer(I4B),
intent(in) :: iout
222 integer(I4B) :: varid
223 logical(LGP) :: layered
225 layered = (idt%layered .and.
is_layered(nc_vars%grid))
231 varid = nc_vars%varid(idt%tagname)
241 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: int1d
243 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
246 integer(I4B),
intent(in) :: varid
247 character(len=*),
intent(in) :: input_fname
248 integer(I4B),
dimension(:),
allocatable :: array_shape
249 integer(I4B),
dimension(:, :, :),
contiguous,
pointer :: int3d_ptr
250 integer(I4B),
dimension(:, :),
contiguous,
pointer :: int2d_ptr
251 integer(I4B) :: nvals
256 if (idt%shape ==
'NODES')
then
258 nvals = product(mshape)
259 if (
size(mshape) == 3)
then
260 int3d_ptr(1:mshape(3), 1:mshape(2), 1:mshape(1)) => int1d(1:nvals)
261 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int3d_ptr), &
263 else if (
size(mshape) == 2)
then
264 int2d_ptr(1:mshape(2), 1:mshape(1)) => int1d(1:nvals)
265 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int2d_ptr), &
267 else if (
size(mshape) == 1)
then
268 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d), nc_vars%nc_fname)
274 nvals = array_shape(1)
276 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d), nc_vars%nc_fname)
286 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: int1d
288 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
291 integer(I4B),
intent(in) :: iper
292 character(len=*),
intent(in) :: input_fname
293 integer(I4B),
dimension(:),
allocatable :: layer_shape
294 integer(I4B) :: varid, nlay, ncpl, istp
299 varid = nc_vars%varid(idt%tagname)
302 ncpl = product(layer_shape)
304 if (
size(mshape) == 3)
then
305 select case (idt%shape)
306 case (
'NCPL',
'NAUX NCPL')
307 if (nc_vars%grid ==
'STRUCTURED')
then
308 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d, &
309 start=(/1, 1, istp/), &
310 count=(/mshape(3), mshape(2), 1/)), &
312 else if (nc_vars%grid ==
'LAYERED MESH')
then
313 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d, &
314 start=(/1, istp/), count=(/ncpl, 1/)), &
317 case (
'NODES',
'NAUX NODES')
318 write (
errmsg,
'(a,a,a)') &
319 'Timeseries netcdf input read not supported for DIS full grid int1d &
320 &type ('//trim(idt%tagname)//
').'
332 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: int1d
334 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
337 character(len=*),
intent(in) :: input_fname
338 integer(I4B),
dimension(:),
allocatable :: layer_shape
339 integer(I4B) :: nlay, varid
340 integer(I4B) :: k, ncpl
341 integer(I4B) :: index_start, index_stop
342 integer(I4B),
dimension(:),
contiguous,
pointer :: int1d_ptr
348 ncpl = product(layer_shape)
351 varid = nc_vars%varid(idt%tagname, layer=k)
352 index_stop = index_start + ncpl - 1
353 int1d_ptr(1:ncpl) => int1d(index_start:index_stop)
354 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
356 index_start = index_stop + 1
366 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: int1d
368 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
371 integer(I4B),
intent(in) :: iper
372 character(len=*),
intent(in) :: input_fname
373 integer(I4B),
dimension(:),
allocatable :: layer_shape
374 integer(I4B) :: nlay, varid
375 integer(I4B) :: ncpl, nvals, istp
380 nvals = product(mshape)
381 ncpl = product(layer_shape)
383 varid = nc_vars%varid(idt%tagname)
384 select case (idt%shape)
385 case (
'NCPL',
'NAUX NCPL')
386 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d, &
387 start=(/1, istp/), count=(/ncpl, 1/)), &
389 case (
'NODES',
'NAUX NODES')
390 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d, &
391 start=(/1, istp/), count=(/nvals, 1/)), &
401 integer(I4B),
dimension(:, :),
contiguous,
pointer,
intent(in) :: int2d
403 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
406 integer(I4B),
intent(in) :: varid
407 character(len=*),
intent(in) :: input_fname
408 integer(I4B),
dimension(:),
contiguous,
pointer :: int1d_ptr
409 integer(I4B),
dimension(:),
allocatable :: array_shape
410 integer(I4B) :: ncpl, nlay
414 if (nc_vars%grid ==
'STRUCTURED')
then
415 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int2d), nc_vars%nc_fname)
416 else if (nc_vars%grid ==
'LAYERED MESH')
then
418 ncpl = product(array_shape)
419 int1d_ptr(1:ncpl) => int2d(:, :)
420 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
429 integer(I4B),
dimension(:, :),
contiguous,
pointer,
intent(in) :: int2d
431 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
434 character(len=*),
intent(in) :: input_fname
435 integer(I4B),
dimension(:),
allocatable :: layer_shape
437 integer(I4B) :: ncpl, nlay, varid
438 integer(I4B),
dimension(:),
contiguous,
pointer :: int1d_ptr
442 if (
size(mshape) == 3)
then
443 write (
errmsg,
'(a,a,a)') &
444 'Layered netcdf read not supported for DIS int2d type ('// &
445 trim(idt%tagname)//
').'
448 else if (
size(mshape) == 2)
then
450 ncpl = layer_shape(1)
452 varid = nc_vars%varid(idt%tagname, layer=k)
453 int1d_ptr(1:ncpl) => int2d(1:ncpl, k)
454 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
464 integer(I4B),
dimension(:, :, :),
contiguous,
pointer,
intent(in) :: int3d
466 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
469 integer(I4B),
intent(in) :: varid
470 character(len=*),
intent(in) :: input_fname
471 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int3d), nc_vars%nc_fname)
478 integer(I4B),
dimension(:, :, :),
contiguous,
pointer,
intent(in) :: int3d
480 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
483 character(len=*),
intent(in) :: input_fname
484 integer(I4B),
dimension(:),
allocatable :: layer_shape
486 integer(I4B) :: ncpl, nlay, varid
487 integer(I4B) :: index_start, index_stop
488 integer(I4B),
dimension(:),
contiguous,
pointer :: int1d_ptr
493 ncpl = product(layer_shape)
496 varid = nc_vars%varid(idt%tagname, layer=k)
497 index_stop = index_start + ncpl - 1
498 int1d_ptr(1:ncpl) => int3d(:, :, k:k)
499 call nf_verify(nf90_get_var(nc_vars%ncid, varid, int1d_ptr), &
501 index_start = index_stop + 1
509 real(DP),
dimension(:),
contiguous,
pointer,
intent(in) :: dbl1d
511 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
514 integer(I4B),
intent(in) :: varid
515 character(len=*),
intent(in) :: input_fname
516 integer(I4B),
dimension(:),
allocatable :: array_shape
517 real(DP),
dimension(:, :, :),
contiguous,
pointer :: dbl3d_ptr
518 real(DP),
dimension(:, :),
contiguous,
pointer :: dbl2d_ptr
519 integer(I4B) :: nvals
524 if (idt%shape ==
'NODES')
then
526 nvals = product(mshape)
527 if (
size(mshape) == 3)
then
528 dbl3d_ptr(1:mshape(3), 1:mshape(2), 1:mshape(1)) => dbl1d(1:nvals)
529 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl3d_ptr), &
531 else if (
size(mshape) == 2)
then
532 dbl2d_ptr(1:mshape(2), 1:mshape(1)) => dbl1d(1:nvals)
533 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl2d_ptr), &
535 else if (
size(mshape) == 1)
then
536 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d), nc_vars%nc_fname)
542 nvals = array_shape(1)
544 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d), nc_vars%nc_fname)
551 iper, input_fname, iaux)
554 real(DP),
dimension(:),
contiguous,
pointer,
intent(in) :: dbl1d
556 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
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 real(DP),
dimension(:, :, :),
contiguous,
pointer :: dbl3d
564 integer(I4B) :: varid, nlay, ncpl, nvals
565 integer(I4B) :: n, istp
572 if (
present(iaux))
then
573 varid = nc_vars%varid(idt%tagname, iaux=iaux)
575 varid = nc_vars%varid(idt%tagname)
579 ncpl = product(layer_shape)
580 nvals = product(mshape)
582 if (
size(mshape) == 3)
then
583 select case (idt%shape)
584 case (
'NCPL',
'NAUX NCPL')
585 if (nc_vars%grid ==
'STRUCTURED')
then
586 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d, &
587 start=(/1, 1, istp/), &
588 count=(/mshape(3), mshape(2), 1/)), &
590 else if (nc_vars%grid ==
'LAYERED MESH')
then
591 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d, &
592 start=(/1, istp/), count=(/ncpl, 1/)), &
595 case (
'NODES',
'NAUX NODES')
596 if (nc_vars%grid ==
'STRUCTURED')
then
597 dbl3d(1:mshape(3), 1:mshape(2), 1:mshape(1)) => dbl1d(1:nvals)
598 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl3d, &
599 start=(/1, 1, 1, istp/), &
600 count=(/mshape(3), mshape(2), mshape(1), &
601 1/)), nc_vars%nc_fname)
602 else if (nc_vars%grid ==
'LAYERED MESH')
then
603 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d, &
604 start=(/1, istp/), count=(/nvals, 1/)), &
616 real(DP),
dimension(:),
contiguous,
pointer,
intent(in) :: dbl1d
618 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
621 character(len=*),
intent(in) :: input_fname
622 integer(I4B),
dimension(:),
allocatable :: layer_shape
623 integer(I4B) :: nlay, varid
624 integer(I4B) :: k, ncpl
625 integer(I4B) :: index_start, index_stop
626 real(DP),
dimension(:),
contiguous,
pointer :: dbl1d_ptr
631 ncpl = product(layer_shape)
634 varid = nc_vars%varid(idt%tagname, layer=k)
635 index_stop = index_start + ncpl - 1
636 dbl1d_ptr(1:ncpl) => dbl1d(index_start:index_stop)
637 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
639 index_start = index_stop + 1
646 iper, input_fname, iaux)
649 real(DP),
dimension(:),
contiguous,
pointer,
intent(in) :: dbl1d
651 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
654 integer(I4B),
intent(in) :: iper
655 character(len=*),
intent(in) :: input_fname
656 integer(I4B),
optional,
intent(in) :: iaux
657 integer(I4B),
dimension(:),
allocatable :: layer_shape
658 integer(I4B) :: nlay, varid
659 integer(I4B) :: k, n, ncpl, idx, istp
660 real(DP),
dimension(:),
contiguous,
pointer :: dbl1d_ptr
665 ncpl = product(layer_shape)
666 allocate (dbl1d_ptr(ncpl))
669 if (
present(iaux))
then
670 varid = nc_vars%varid(idt%tagname, layer=k, iaux=iaux)
672 varid = nc_vars%varid(idt%tagname, layer=k)
674 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr, &
675 start=(/1, istp/), count=(/ncpl, 1/)), &
677 if (idt%shape ==
'NODES' .or. idt%shape ==
'NAUX NODES')
then
679 idx = (k - 1) * ncpl + n
680 dbl1d(idx) = dbl1d_ptr(n)
682 else if (idt%shape ==
'NCPL' .or. idt%shape ==
'NAUX NCPL')
then
684 dbl1d(n) = dbl1d_ptr(n)
690 deallocate (dbl1d_ptr)
697 real(DP),
dimension(:, :),
contiguous,
pointer,
intent(in) :: dbl2d
699 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
702 integer(I4B),
intent(in) :: varid
703 character(len=*),
intent(in) :: input_fname
704 real(DP),
dimension(:),
contiguous,
pointer :: dbl1d_ptr
705 integer(I4B),
dimension(:),
allocatable :: array_shape
706 integer(I4B) :: ncpl, nlay
710 if (nc_vars%grid ==
'STRUCTURED')
then
711 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl2d), nc_vars%nc_fname)
712 else if (nc_vars%grid ==
'LAYERED MESH')
then
714 ncpl = product(array_shape)
715 dbl1d_ptr(1:ncpl) => dbl2d(:, :)
716 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
725 real(DP),
dimension(:, :),
contiguous,
pointer,
intent(in) :: dbl2d
727 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
730 character(len=*),
intent(in) :: input_fname
731 integer(I4B),
dimension(:),
allocatable :: layer_shape
733 integer(I4B) :: ncpl, nlay, varid
734 real(DP),
dimension(:),
contiguous,
pointer :: dbl1d_ptr
738 if (
size(mshape) == 3)
then
739 write (
errmsg,
'(a,a,a)') &
740 'Layered netcdf read not supported for DIS dbl2d type ('// &
741 trim(idt%tagname)//
').'
744 else if (
size(mshape) == 2)
then
746 ncpl = layer_shape(1)
748 varid = nc_vars%varid(idt%tagname, layer=k)
749 dbl1d_ptr(1:ncpl) => dbl2d(1:ncpl, k)
750 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
760 real(DP),
dimension(:, :, :),
contiguous,
pointer,
intent(in) :: dbl3d
762 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
765 integer(I4B),
intent(in) :: varid
766 character(len=*),
intent(in) :: input_fname
768 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl3d), nc_vars%nc_fname)
775 real(DP),
dimension(:, :, :),
contiguous,
pointer,
intent(in) :: dbl3d
777 integer(I4B),
dimension(:),
contiguous,
pointer,
intent(in) :: mshape
780 character(len=*),
intent(in) :: input_fname
781 integer(I4B),
dimension(:),
allocatable :: layer_shape
783 integer(I4B) :: ncpl, nlay, varid
784 integer(I4B) :: index_start, index_stop
785 real(DP),
dimension(:),
contiguous,
pointer :: dbl1d_ptr
791 ncpl = product(layer_shape)
794 varid = nc_vars%varid(idt%tagname, layer=k)
795 index_stop = index_start + ncpl - 1
796 dbl1d_ptr(1:ncpl) => dbl3d(:, :, k:k)
797 call nf_verify(nf90_get_var(nc_vars%ncid, varid, dbl1d_ptr), &
799 index_start = index_stop + 1
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
real(dp), parameter dnodata
real no data constant
This module defines variable data types.
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_integer1d_layered_spd(int1d, mf6_input, mshape, idt, nc_vars, iper, input_fname)
load type 1d integer layered
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_integer1d_spd(int1d, mf6_input, mshape, idt, nc_vars, iper, input_fname)
load type 1d 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.
This module contains the NetCDFCommonModule.
integer(i4b) function, public ixstp()
step index for timeseries data
subroutine, public nf_verify(res, nc_fname)
error check a netcdf-fortran interface call
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
This module contains the SourceCommonModule.
subroutine, public get_layered_shape(mshape, nlay, layer_shape)
subroutine, public get_shape_from_string(shape_string, array_shape, memoryPath)
Type describing input variables for a package in NetCDF file.