31 character(len=LENMEMPATH) :: memorypath
32 character(len=LENMEMPATH) :: input_mempath =
''
33 character(len=LENMODELNAME),
pointer :: name_model => null()
34 character(len=LINELENGTH),
pointer :: input_fname => null()
35 character(len=LINELENGTH),
pointer :: output_fname => null()
36 integer(I4B),
pointer :: inunit => null()
37 integer(I4B),
pointer :: iout => null()
38 integer(I4B),
pointer :: nodes => null()
39 integer(I4B),
pointer :: nodesuser => null()
40 integer(I4B),
pointer :: nja => null()
41 integer(I4B),
pointer :: njas => null()
42 integer(I4B),
pointer :: lenuni => null()
43 integer(I4B),
pointer :: ndim => null()
44 integer(I4B),
pointer :: icondir => null()
45 integer(I4B),
pointer :: nogrb => null()
46 real(dp),
dimension(:),
pointer,
contiguous :: xc => null()
47 real(dp),
dimension(:),
pointer,
contiguous :: yc => null()
48 real(dp),
pointer :: yorigin => null()
49 real(dp),
pointer :: xorigin => null()
50 real(dp),
pointer :: angrot => null()
51 integer(I4B),
dimension(:),
pointer,
contiguous :: mshape => null()
52 real(dp),
dimension(:),
pointer,
contiguous :: top => null()
53 real(dp),
dimension(:),
pointer,
contiguous :: bot => null()
54 real(dp),
dimension(:),
pointer,
contiguous :: area => null()
57 real(dp),
dimension(:),
pointer,
contiguous :: dbuff => null()
58 integer(I4B),
dimension(:),
pointer,
contiguous :: ibuff => null()
59 integer(I4B),
dimension(:),
pointer,
contiguous :: nodereduced => null()
60 integer(I4B),
dimension(:),
pointer,
contiguous :: nodeuser => null()
130 call store_error(
'Programmer error: dis_df must be overridden', &
140 integer(I4B),
intent(in) :: moffset
143 integer(I4B) :: i, j, ipos, iglo, jglo
146 do ipos = this%con%ia(i), this%con%ia(i + 1) - 1
147 j = this%con%ja(ipos)
150 call sparse%addconnection(iglo, jglo, 1)
156 subroutine dis_mc(this, moffset, idxglo, matrix_sln)
159 integer(I4B),
intent(in) :: moffset
160 integer(I4B),
dimension(:),
intent(inout) :: idxglo
163 integer(I4B) :: i, j, ipos, iglo, jglo
167 do ipos = this%con%ia(i), this%con%ia(i + 1) - 1
168 j = this%con%ja(ipos)
170 idxglo(ipos) = matrix_sln%get_position(iglo, jglo)
179 integer(I4B),
dimension(:),
intent(in) :: icelltype
181 integer(I4B),
dimension(:),
allocatable :: ict
182 integer(I4B) :: nu, nr
185 allocate (ict(this%nodesuser))
186 do nu = 1, this%nodesuser
187 nr = this%get_nodenumber(nu, 0)
189 ict(nu) = icelltype(nr)
195 if (this%nogrb == 0)
call this%write_grb(ict)
201 integer(I4B),
dimension(:),
intent(in) :: icelltype
202 call store_error(
'Programmer error: write_grb must be overridden', &
214 deallocate (this%name_model)
215 deallocate (this%input_fname)
216 deallocate (this%output_fname)
244 call this%con%con_da()
245 deallocate (this%con)
251 integer(I4B),
intent(in) :: nodeu
252 character(len=*),
intent(inout) :: str
254 call store_error(
'Programmer error: nodeu_to_string must be overridden', &
261 integer(I4B),
intent(in) :: nodeu
262 integer(I4B),
dimension(:),
intent(inout) :: arr
264 call store_error(
'Programmer error: nodeu_to_array must be overridden', &
271 integer(I4B),
intent(in) :: noder
272 integer(I4B) :: nodenumber
274 if (this%nodes < this%nodesuser)
then
275 nodenumber = this%nodeuser(noder)
283 integer(I4B),
intent(in) :: nodeu
284 integer(I4B),
intent(in) :: icheck
285 integer(I4B) :: nodenumber
288 call store_error(
'Programmer error: get_nodenumber_idx1 must be overridden', &
294 integer(I4B),
intent(in) :: k, j
295 integer(I4B),
intent(in) :: icheck
296 integer(I4B) :: nodenumber
299 call store_error(
'Programmer error: get_nodenumber_idx2 must be overridden', &
305 integer(I4B),
intent(in) :: k, i, j
306 integer(I4B),
intent(in) :: icheck
307 integer(I4B) :: nodenumber
310 call store_error(
'Programmer error: get_nodenumber_idx3 must be overridden', &
319 integer(I4B),
intent(in) :: noden
320 integer(I4B),
intent(in) :: nodem
321 integer(I4B),
intent(in) :: ihc
322 real(DP),
intent(inout) :: xcomp
323 real(DP),
intent(inout) :: ycomp
324 real(DP),
intent(inout) :: zcomp
325 integer(I4B),
intent(in) :: ipos
327 call store_error(
'Programmer error: connection_normal must be overridden', &
335 xcomp, ycomp, zcomp, conlen)
337 integer(I4B),
intent(in) :: noden
338 integer(I4B),
intent(in) :: nodem
339 logical,
intent(in) :: nozee
340 real(DP),
intent(in) :: satn
341 real(DP),
intent(in) :: satm
342 integer(I4B),
intent(in) :: ihc
343 real(DP),
intent(inout) :: xcomp
344 real(DP),
intent(inout) :: ycomp
345 real(DP),
intent(inout) :: zcomp
346 real(DP),
intent(inout) :: conlen
348 call store_error(
'Programmer error: connection_vector must be overridden', &
354 real(dp),
intent(in) :: x
355 real(dp),
intent(in) :: y
356 real(dp),
intent(in) :: xorigin
357 real(dp),
intent(in) :: yorigin
358 real(dp),
intent(in) :: angrot
359 real(dp),
intent(out) :: xglo
360 real(dp),
intent(out) :: yglo
369 if (ang /=
dzero)
then
370 xglo = x * cos(ang) - y * sin(ang)
371 yglo = x * sin(ang) + y * cos(ang)
375 xglo = xglo + xorigin
376 yglo = yglo + yorigin
382 character(len=*),
intent(out) :: dis_type
384 dis_type =
"Not implemented"
385 call store_error(
'Programmer error: get_dis_type must be overridden', &
393 integer(I4B) :: dis_enum
396 call store_error(
'Programmer error: get_dis_enum must be overridden', &
404 character(len=*),
intent(in) :: name_model
405 character(len=*),
intent(in) :: input_mempath
406 logical(LGP) :: found
412 allocate (this%name_model)
413 allocate (this%input_fname)
414 allocate (this%output_fname)
416 call mem_allocate(this%inunit,
'INUNIT', this%memoryPath)
419 call mem_allocate(this%nodesuser,
'NODESUSER', this%memoryPath)
421 call mem_allocate(this%icondir,
'ICONDIR', this%memoryPath)
423 call mem_allocate(this%xorigin,
'XORIGIN', this%memoryPath)
424 call mem_allocate(this%yorigin,
'YORIGIN', this%memoryPath)
425 call mem_allocate(this%angrot,
'ANGROT', this%memoryPath)
428 call mem_allocate(this%lenuni,
'LENUNI', this%memoryPath)
431 this%name_model = name_model
432 this%input_mempath = input_mempath
433 this%input_fname =
''
434 this%output_fname =
''
451 this%input_mempath, found)
453 this%input_mempath, found)
454 if (.not. found)
then
455 this%output_fname = trim(this%input_fname)//
'.grb'
465 call mem_allocate(this%mshape, this%ndim,
'MSHAPE', this%memoryPath)
466 call mem_allocate(this%xc, this%nodes,
'XC', this%memoryPath)
467 call mem_allocate(this%yc, this%nodes,
'YC', this%memoryPath)
468 call mem_allocate(this%top, this%nodes,
'TOP', this%memoryPath)
469 call mem_allocate(this%bot, this%nodes,
'BOT', this%memoryPath)
470 call mem_allocate(this%area, this%nodes,
'AREA', this%memoryPath)
473 this%mshape(1) = this%nodes
476 if (this%nodes < this%nodesuser)
then
477 isize = this%nodesuser
483 call mem_allocate(this%dbuff, isize,
'DBUFF', this%name_model)
484 call mem_allocate(this%ibuff, isize,
'IBUFF', this%name_model)
494 flag_string, allow_zero)
result(nodeu)
497 integer(I4B),
intent(inout) :: lloc
498 integer(I4B),
intent(inout) :: istart
499 integer(I4B),
intent(inout) :: istop
500 integer(I4B),
intent(in) :: in
501 integer(I4B),
intent(in) :: iout
502 character(len=*),
intent(inout) :: line
503 logical,
optional,
intent(in) :: flag_string
504 logical,
optional,
intent(in) :: allow_zero
505 integer(I4B) :: nodeu
508 call store_error(
'Programmer error: nodeu_from_string must be overridden', &
521 allow_zero)
result(nodeu)
524 character(len=*),
intent(inout) :: cellid
525 integer(I4B),
intent(in) :: inunit
526 integer(I4B),
intent(in) :: iout
527 logical,
optional,
intent(in) :: flag_string
528 logical,
optional,
intent(in) :: allow_zero
529 integer(I4B) :: nodeu
532 call store_error(
'Programmer error: nodeu_from_cellid must be overridden', &
544 flag_string)
result(noder)
547 integer(I4B),
intent(inout) :: lloc
548 integer(I4B),
intent(inout) :: istart
549 integer(I4B),
intent(inout) :: istop
550 integer(I4B),
intent(in) :: in
551 integer(I4B),
intent(in) :: iout
552 character(len=*),
intent(inout) :: line
553 logical,
optional,
intent(in) :: flag_string
554 integer(I4B) :: noder
556 integer(I4B) :: nodeu
557 character(len=LINELENGTH) :: nodestr
558 logical :: flag_string_local
560 if (
present(flag_string))
then
561 flag_string_local = flag_string
563 flag_string_local = .false.
565 nodeu = this%nodeu_from_string(lloc, istart, istop, in, iout, line, &
570 noder = this%get_nodenumber(nodeu, 0)
574 if (noder <= 0 .and. .not. flag_string_local)
then
575 call this%nodeu_to_string(nodeu, nodestr)
577 ' Cell is outside active grid domain: '// &
578 trim(adjustl(nodestr))
592 allow_zero)
result(noder)
594 integer(I4B) :: noder
597 character(len=*),
intent(inout) :: cellid
598 integer(I4B),
intent(in) :: inunit
599 integer(I4B),
intent(in) :: iout
600 logical,
optional,
intent(in) :: flag_string
601 logical,
optional,
intent(in) :: allow_zero
603 integer(I4B) :: nodeu
604 logical :: allowzerolocal
605 character(len=LINELENGTH) :: nodestr
606 logical :: flag_string_local
608 if (
present(flag_string))
then
609 flag_string_local = flag_string
611 flag_string_local = .false.
613 if (
present(allow_zero))
then
614 allowzerolocal = allow_zero
616 allowzerolocal = .false.
619 nodeu = this%nodeu_from_cellid(cellid, inunit, iout, flag_string_local, &
624 noder = this%get_nodenumber(nodeu, 0)
628 if (noder <= 0 .and. .not. flag_string_local)
then
629 call this%nodeu_to_string(nodeu, nodestr)
631 ' Cell is outside active grid domain: '// &
632 trim(adjustl(nodestr))
641 call store_error(
'Programmer error: supports_layers must be overridden', &
651 call store_error(
'Programmer error: get_ncpl must be overridden', &
661 integer(I4B),
intent(in) :: n
662 real(dp),
intent(in) :: x
673 thick = (tp - bt) * sat
681 integer(I4B),
intent(in) :: ic
682 real(DP),
allocatable,
intent(out) :: polyverts(:, :)
683 logical(LGP),
intent(in),
optional :: closed
685 errmsg =
'Programmer error: get_polyverts must be overridden'
694 character(len=*),
intent(inout) :: line
695 integer(I4B),
intent(inout) :: lloc
696 integer(I4B),
intent(inout) :: istart
697 integer(I4B),
intent(inout) :: istop
698 integer(I4B),
intent(in) :: in
699 integer(I4B),
intent(in) :: iout
700 integer(I4B),
dimension(:),
pointer,
contiguous,
intent(inout) :: iarray
701 character(len=*),
intent(in) :: aname
703 errmsg =
'Programmer error: read_int_array must be overridden'
712 character(len=*),
intent(inout) :: line
713 integer(I4B),
intent(inout) :: lloc
714 integer(I4B),
intent(inout) :: istart
715 integer(I4B),
intent(inout) :: istop
716 integer(I4B),
intent(in) :: in
717 integer(I4B),
intent(in) :: iout
718 real(DP),
dimension(:),
pointer,
contiguous,
intent(inout) :: darray
719 character(len=*),
intent(in) :: aname
721 errmsg =
'Programmer error: read_dbl_array must be overridden'
729 integer(I4B),
dimension(:),
pointer,
contiguous,
intent(in) :: ibuff1
730 integer(I4B),
dimension(:),
pointer,
contiguous,
intent(inout) :: ibuff2
732 integer(I4B) :: nodeu
733 integer(I4B) :: noder
735 do nodeu = 1, this%nodesuser
736 noder = this%get_nodenumber(nodeu, 0)
737 if (noder <= 0) cycle
738 ibuff2(noder) = ibuff1(nodeu)
746 real(DP),
dimension(:),
pointer,
contiguous,
intent(in) :: buff1
747 real(DP),
dimension(:),
pointer,
contiguous,
intent(inout) :: buff2
749 integer(I4B) :: nodeu
750 integer(I4B) :: noder
752 do nodeu = 1, this%nodesuser
753 noder = this%get_nodenumber(nodeu, 0)
754 if (noder <= 0) cycle
755 buff2(noder) = buff1(nodeu)
766 subroutine read_list(this, line_reader, in, iout, iprpak, nlist, &
767 inamedbound, iauxmultcol, nodelist, rlist, auxvar, &
768 auxname, boundname, label, pkgname, tsManager, iscloc, &
781 integer(I4B),
intent(in) :: in
782 integer(I4B),
intent(in) :: iout
783 integer(I4B),
intent(in) :: iprpak
784 integer(I4B),
intent(inout) :: nlist
785 integer(I4B),
intent(in) :: inamedbound
786 integer(I4B),
intent(in) :: iauxmultcol
787 integer(I4B),
dimension(:),
pointer,
contiguous,
intent(inout) :: nodelist
788 real(DP),
dimension(:, :),
pointer,
contiguous,
intent(inout) :: rlist
789 real(DP),
dimension(:, :),
pointer,
contiguous,
intent(inout) :: auxvar
790 character(len=LENAUXNAME),
dimension(:),
intent(inout) :: auxname
791 character(len=LENBOUNDNAME),
dimension(:),
pointer,
contiguous, &
792 intent(inout) :: boundname
793 character(len=*),
intent(in) :: label
794 character(len=*),
intent(in) :: pkgName
796 integer(I4B),
intent(in) :: iscloc
797 integer(I4B),
intent(in),
optional :: indxconvertflux
800 integer(I4B) :: nodeu, noder
801 character(len=LINELENGTH) :: nodestr
802 integer(I4B) :: ii, jj
803 real(DP),
pointer :: bndElem => null()
809 call lstrdobj%read_list(line_reader, in, iout, nlist, inamedbound, &
810 this%mshape, nodelist, rlist, auxvar, auxname, &
815 if (lstrdobj%ntxtrlist > 0)
then
816 do l = 1, lstrdobj%ntxtrlist
817 ii = lstrdobj%idxtxtrow(l)
818 jj = lstrdobj%idxtxtcol(l)
820 bndelem => rlist(jj, ii)
822 pkgname,
'BND', tsmanager, iprpak, &
824 if (
associated(tslinkbnd))
then
829 if (iauxmultcol > 0 .and. jj == iscloc)
then
830 tslinkbnd%RMultiplier => auxvar(iauxmultcol, ii)
834 if (lstrdobj%inamedbound == 1)
then
835 tslinkbnd%BndName = lstrdobj%boundname(tslinkbnd%IRow)
840 if (
present(indxconvertflux))
then
841 if (indxconvertflux == jj)
then
842 tslinkbnd%convertflux = .true.
844 noder = this%get_nodenumber(nodeu, 0)
845 tslinkbnd%CellArea = this%get_area(noder)
854 if (lstrdobj%ntxtauxvar > 0)
then
855 do l = 1, lstrdobj%ntxtauxvar
856 ii = lstrdobj%idxtxtauxrow(l)
857 jj = lstrdobj%idxtxtauxcol(l)
859 bndelem => auxvar(jj, ii)
861 pkgname,
'AUX', tsmanager, iprpak, &
863 if (lstrdobj%inamedbound == 1)
then
864 if (
associated(tslinkaux))
then
865 tslinkaux%BndName = lstrdobj%boundname(tslinkaux%IRow)
872 if (iauxmultcol > 0)
then
874 rlist(iscloc, l) = rlist(iscloc, l) * auxvar(iauxmultcol, l)
879 if (iprpak /= 0)
then
880 call lstrdobj%write_list()
886 if (this%nodes < this%nodesuser)
then
889 noder = this%get_nodenumber(nodeu, 0)
891 call this%nodeu_to_string(nodeu, nodestr)
893 ' Cell is outside active grid domain: '// &
894 trim(adjustl(nodestr))
915 icolbnd, aname, inunit, iout)
918 integer(I4B),
intent(in) :: ncolbnd
919 integer(I4B),
intent(in) :: maxbnd
920 integer(I4B),
dimension(maxbnd) :: nodelist
921 real(DP),
dimension(ncolbnd, maxbnd),
intent(inout) :: darray
922 integer(I4B),
intent(in) :: icolbnd
923 character(len=*),
intent(in) :: aname
924 integer(I4B),
intent(in) :: inunit
925 integer(I4B),
intent(in) :: iout
927 errmsg =
'Programmer error: read_layer_array must be overridden'
936 cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
939 real(DP),
dimension(:),
pointer,
contiguous,
intent(inout) :: darray
940 integer(I4B),
intent(in) :: iout
941 integer(I4B),
intent(in) :: iprint
942 integer(I4B),
intent(in) :: idataun
943 character(len=*),
intent(in) :: aname
944 character(len=*),
intent(in) :: cdatafmp
945 integer(I4B),
intent(in) :: nvaluesp
946 integer(I4B),
intent(in) :: nwidthp
947 character(len=*),
intent(in) :: editdesc
948 real(DP),
intent(in) :: dinact
950 errmsg =
'Programmer error: record_array must be overridden'
958 real(DP),
dimension(:),
intent(in) :: flowja
959 integer(I4B),
intent(in) :: ibinun
960 integer(I4B),
intent(in) :: iout
962 character(len=16),
dimension(1) :: text
964 data text(1)/
' FLOW-JA-FACE'/
967 call ubdsv1(
kstp,
kper, text(1), ibinun, flowja,
size(flowja), 1, 1, &
975 integer(I4B),
intent(in) :: noder
976 character(len=*),
intent(inout) :: str
978 integer(I4B) :: nodeu
980 nodeu = this%get_nodeuser(noder)
981 call this%nodeu_to_string(nodeu, str)
988 integer(I4B),
intent(in) :: noder
989 integer(I4B),
dimension(:),
intent(inout) :: arr
991 integer(I4B) :: nodeu
993 nodeu = this%get_nodeuser(noder)
994 call this%nodeu_to_array(nodeu, arr)
999 dstmodel, dstpackage, naux, auxtxt, &
1000 ibdchn, nlist, iout)
1002 character(len=16),
intent(in) :: text
1003 character(len=16),
intent(in) :: textmodel
1004 character(len=16),
intent(in) :: textpackage
1005 character(len=16),
intent(in) :: dstmodel
1006 character(len=16),
intent(in) :: dstpackage
1007 integer(I4B),
intent(in) :: naux
1008 character(len=16),
dimension(:),
intent(in) :: auxtxt
1009 integer(I4B),
intent(in) :: ibdchn
1010 integer(I4B),
intent(in) :: nlist
1011 integer(I4B),
intent(in) :: iout
1013 errmsg =
'Programmer error: record_srcdst_list_header must be overridden'
1019 naux, aux, olconv, olconv2)
1022 integer(I4B),
intent(in) :: ibdchn
1023 integer(I4B),
intent(in) :: noder
1024 integer(I4B),
intent(in) :: noder2
1025 real(DP),
intent(in) :: q
1026 integer(I4B),
intent(in) :: naux
1027 real(DP),
dimension(naux),
intent(in) :: aux
1028 logical,
optional,
intent(in) :: olconv
1029 logical,
optional,
intent(in) :: olconv2
1033 integer(I4B) :: nodeu
1034 integer(I4B) :: nodeu2
1037 if (
present(olconv))
then
1043 nodeu = this%get_nodeuser(noder)
1047 if (
present(olconv2))
then
1053 nodeu2 = this%get_nodeuser(noder2)
1057 call ubdsvd(ibdchn, nodeu, nodeu2, q, naux, aux)
1066 integer(I4B),
intent(in) :: maxbnd
1067 integer(I4B),
dimension(:),
pointer,
contiguous :: darray
1068 integer(I4B),
dimension(maxbnd),
intent(inout) :: nodelist
1069 integer(I4B),
intent(inout) :: nbound
1070 character(len=*),
intent(in) :: aname
1072 errmsg =
'Programmer error: nlarray_to_nodelist must be overridden'
1080 integer(I4B),
intent(inout) :: n
1081 integer(I4B),
dimension(:),
intent(in) :: ibound
1083 integer(I4B) :: m, ii, iis
1084 logical done, bottomcell
1089 do while (.not. done)
1091 cloop:
do ii = this%con%ia(n) + 1, this%con%ia(n + 1) - 1
1093 iis = this%con%jas(ii)
1094 if (this%con%ihc(iis) == 0 .and. m > n)
then
1097 bottomcell = .false.
1100 if (ibound(m) /= 0)
then
1110 if (bottomcell) done = .true.
1117 integer(I4B),
intent(in) :: node
1120 area = this%area(node)
1134 real(dp) :: area_factor
1137 integer(I4B),
intent(in) :: node
1138 integer(I4B),
intent(in) :: idx_conn
1140 real(dp) :: area_node
1141 real(dp) :: area_conn
1144 area_node = this%area(node)
1145 area_conn = this%con%hwva(idx_conn)
1148 area_factor = area_conn / area_node
1161 integer(I4B),
intent(in) :: n
1162 integer(I4B),
intent(in) :: m
1163 integer(I4B),
intent(in) :: idx_conn
1164 real(DP),
intent(out) :: width_n
1165 real(DP),
intent(out) :: width_m
1167 integer(I4B) :: isympos
1170 isympos = this%con%jas(idx_conn)
1171 width_n = this%con%hwva(isympos)
1183 select case (this%get_dis_enum())
1196 select case (this%get_dis_enum())
1209 select case (this%get_dis_enum())
real(dp) function get_area(this, node)
Return the cell area for the given node.
subroutine dis_df(this)
Define the discretization.
subroutine noder_to_string(this, noder, str)
Convert reduced node number to string (nodenumber), (k,j) or (k,i,j)
integer(i4b) function get_nodenumber_idx2(this, k, j, icheck)
real(dp) function get_area_factor(this, node, idx_conn)
@ brief Calculate the area factor for the cell connection
subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, icolbnd, aname, inunit, iout)
Read a 2d double array into col icolbnd of darray.
integer(i4b) function nodeu_from_string(this, lloc, istart, istop, in, iout, line, flag_string, allow_zero)
Convert a string to a user nodenumber.
subroutine allocate_scalars(this, name_model, input_mempath)
Allocate and initialize scalar variables.
subroutine dis_mc(this, moffset, idxglo, matrix_sln)
Map cell connections in the numerical solution coefficient matrix.
logical(lgp) function is_1d(this)
@Brief return true if grid is one dimensional
subroutine highest_active(this, n, ibound)
Find the first highest active cell beneath cell n.
subroutine record_array(this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
Record a double precision array.
subroutine record_srcdst_list_header(this, text, textmodel, textpackage, dstmodel, dstpackage, naux, auxtxt, ibdchn, nlist, iout)
Record list header for imeth=6.
integer(i4b) function get_nodenumber_idx3(this, k, i, j, icheck)
subroutine record_srcdst_list_entry(this, ibdchn, noder, noder2, q, naux, aux, olconv, olconv2)
Record list header.
subroutine write_grb(this, icelltype)
Write a binary grid file.
subroutine, public dis_transform_xy(x, y, xorigin, yorigin, angrot, xglo, yglo)
Get global (x, y) coordinates from cell-local coordinates.
subroutine read_int_array(this, line, lloc, istart, istop, iout, in, iarray, aname)
Read an integer array.
logical(lgp) function is_3d(this)
@Brief return true if grid is three dimensional
logical function supports_layers(this)
Indicates whether the grid discretization supports layers.
subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, ipos)
Get normal vector components between the cell and a given neighbor. The normal points outward from th...
subroutine nodeu_to_string(this, nodeu, str)
Convert a user nodenumber to a string (nodenumber), (k,j), or (k,i,j)
integer(i4b) function noder_from_string(this, lloc, istart, istop, in, iout, line, flag_string)
Convert a string to a reduced nodenumber.
subroutine fill_int_array(this, ibuff1, ibuff2)
Fill an integer array.
subroutine get_polyverts(this, ic, polyverts, closed)
Get a 2D array of polygon vertices, listed in clockwise order beginning with the lower left corner.
subroutine dis_da(this)
@brier Deallocate variables
subroutine fill_dbl_array(this, buff1, buff2)
Fill a double precision array.
subroutine dis_ac(this, moffset, sparse)
Add connections to sparse cell connectivity matrix.
subroutine record_connection_array(this, flowja, ibinun, iout)
Record a connection-based double precision array.
subroutine nodeu_to_array(this, nodeu, arr)
Convert a user nodenumber to an array (nodenumber), (k,j), or (k,i,j)
integer(i4b) function get_nodenumber_idx1(this, nodeu, icheck)
integer(i4b) function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, allow_zero)
Convert a cellid string to a user nodenumber.
logical(lgp) function is_2d(this)
@Brief return true if grid is two dimensional
subroutine noder_to_array(this, noder, arr)
Convert reduced node number to array (nodenumber), (k,j) or (k,i,j)
subroutine get_dis_type(this, dis_type)
Get the discretization type (DIS, DISV, or DISU)
integer(i4b) function get_nodeuser(this, noder)
Convert a reduced nodenumber to a user node number.
subroutine get_flow_width(this, n, m, idx_conn, width_n, width_m)
@ brief Calculate the flow width between two cells
subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, darray, aname)
Read a double precision array.
subroutine dis_ar(this, icelltype)
Allocate and setup variables, and write binary grid file.
subroutine allocate_arrays(this)
Allocate and initialize arrays.
integer(i4b) function get_dis_enum(this)
Get the discretization type enumeration.
integer(i4b) function noder_from_cellid(this, cellid, inunit, iout, flag_string, allow_zero)
Convert cellid string to reduced nodenumber.
subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
Convert an integer array to nodelist.
integer(i4b) function get_ncpl(this)
Return number of cells per layer. This is nodes for a DISU grid, as there are no layers.
real(dp) function get_cell_volume(this, n, x)
Return volume of cell n based on x value passed.
subroutine read_list(this, line_reader, in, iout, iprpak, nlist, inamedbound, iauxmultcol, nodelist, rlist, auxvar, auxname, boundname, label, pkgname, tsManager, iscloc, indxconvertflux)
Read a list using the list reader.
subroutine connection_vector(this, noden, nodem, nozee, satn, satm, ihc, xcomp, ycomp, zcomp, conlen)
Get unit vector components between the cell and a given neighbor. Saturation must be provided to comp...
This module contains block parser methods.
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
integer(i4b), parameter lenmodelname
maximum length of the model name
@ disu
DISV6 discretization.
@ dis
DIS6 discretization.
@ dis1d
DIS1D6 discretization.
@ disv2d
DISV2D6 discretization.
@ disv1d
DISV1D6 discretization.
@ dis2d
DIS2D6 discretization.
@ disv
DISU6 discretization.
@ disu1d
DISU1D6 discretization.
@ disundef
undefined discretization
@ disu2d
DISU2D6 discretization.
integer(i4b), parameter lenauxname
maximum length of a aux variable
integer(i4b), parameter lenboundname
maximum length of a bound name
real(dp), parameter dpio180
real constant
real(dp), parameter dzero
real constant zero
integer(i4b), parameter lenmempath
maximum length of the memory path
This module defines variable data types.
Generic List Reader Module.
This module contains the LongLineReaderType.
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
integer(i4b) function, public count_errors()
Return number of errors.
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
real(dp) function squadraticsaturation(top, bot, x, eps)
@ brief sQuadraticSaturation
real(dp), pointer, public pertim
time relative to start of stress period
real(dp), pointer, public totim
time relative to start of simulation
integer(i4b), pointer, public kstp
current time step number
integer(i4b), pointer, public kper
current stress period number
real(dp), pointer, public delt
length of the current time step
subroutine, public read_value_or_time_series(textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, tsLink)
Call this subroutine if the time-series link is available or needed.