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

Data Types

type  disv2dtype
 Vertex grid discretization. More...
 
type  disvfoundtype
 

Functions/Subroutines

subroutine, public disv2d_cr (dis, name_model, input_mempath, inunit, iout)
 Create a new discretization by vertices object. More...
 
subroutine disv2d_load (this)
 Transfer IDM data into this discretization object. More...
 
subroutine disv2d_df (this)
 Define the discretization. More...
 
subroutine disv2d_da (this)
 
subroutine source_options (this)
 Copy options from IDM into package. More...
 
subroutine log_options (this, found)
 Write user options to list file. More...
 
subroutine source_dimensions (this)
 Copy dimensions from IDM into package. More...
 
subroutine log_dimensions (this, found)
 Write dimensions to list file. More...
 
subroutine source_griddata (this)
 Copy grid data from IDM into package. More...
 
subroutine log_griddata (this, found)
 Write griddata found to list file. More...
 
subroutine grid_finalize (this)
 Finalize grid (check properties, allocate arrays, compute connections) More...
 
subroutine source_vertices (this)
 Load grid vertices from IDM into package. More...
 
subroutine define_cellverts (this, icell2d, ncvert, icvert)
 Build data structures to hold cell vertex info. More...
 
subroutine source_cell2d (this)
 Copy cell2d data from IDM into package. More...
 
subroutine connect (this)
 Build grid connections. More...
 
subroutine write_grb (this, icelltype)
 Write a binary grid file. More...
 
subroutine nodeu_to_string (this, nodeu, str)
 Convert a user nodenumber to a string (nodenumber) or (k,j) More...
 
subroutine nodeu_to_array (this, nodeu, arr)
 Convert a user nodenumber to an array (nodenumber) or (k,j) More...
 
integer(i4b) function get_nodenumber_idx1 (this, nodeu, icheck)
 Get reduced node number from user node number. More...
 
subroutine connection_normal (this, noden, nodem, ihc, xcomp, ycomp, zcomp, ipos)
 Get normal vector components between the cell and a given neighbor. More...
 
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. More...
 
subroutine get_dis_type (this, dis_type)
 Get the discretization type. More...
 
integer(i4b) function get_dis_enum (this)
 Get the discretization type enumeration. More...
 
subroutine allocate_scalars (this, name_model, input_mempath)
 Allocate and initialize scalars. More...
 
subroutine allocate_arrays (this)
 Allocate and initialize arrays. More...
 
real(dp) function get_cell2d_area (this, icell2d)
 Get the signed area of the cell. More...
 
integer(i4b) function nodeu_from_string (this, lloc, istart, istop, in, iout, line, flag_string, allow_zero)
 Convert a string to a user nodenumber. More...
 
integer(i4b) function nodeu_from_cellid (this, cellid, inunit, iout, flag_string, allow_zero)
 Convert a cellid string to a user nodenumber. More...
 
subroutine get_polyverts (this, ic, polyverts, closed)
 Get a 2D array of polygon vertices, listed in clockwise order beginning with the lower left corner. More...
 
subroutine record_array (this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
 Record a double precision array. More...
 
subroutine record_srcdst_list_header (this, text, textmodel, textpackage, dstmodel, dstpackage, naux, auxtxt, ibdchn, nlist, iout)
 Record list header for imeth=6. More...
 

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine disv2dmodule::allocate_arrays ( class(disv2dtype this)
private

Definition at line 1018 of file Disv2d.f90.

1019  ! dummy
1020  class(Disv2dType) :: this
1021 
1022  ! Allocate arrays in DisBaseType (mshape, top, bot, area)
1023  call this%DisBaseType%allocate_arrays()
1024  !
1025  ! Allocate arrays for DisvType
1026  if (this%nodes < this%nodesuser) then
1027  call mem_allocate(this%nodeuser, this%nodes, 'NODEUSER', this%memoryPath)
1028  call mem_allocate(this%nodereduced, this%nodesuser, 'NODEREDUCED', &
1029  this%memoryPath)
1030  else
1031  call mem_allocate(this%nodeuser, 1, 'NODEUSER', this%memoryPath)
1032  call mem_allocate(this%nodereduced, 1, 'NODEREDUCED', this%memoryPath)
1033  end if
1034 
1035  ! Initialize
1036  this%mshape(1) = this%nodesuser
1037 

◆ allocate_scalars()

subroutine disv2dmodule::allocate_scalars ( class(disv2dtype this,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  input_mempath 
)

Definition at line 998 of file Disv2d.f90.

999  ! -- dummy
1000  class(Disv2dType) :: this
1001  character(len=*), intent(in) :: name_model
1002  character(len=*), intent(in) :: input_mempath
1003  !
1004  ! -- Allocate parent scalars
1005  call this%DisBaseType%allocate_scalars(name_model, input_mempath)
1006  !
1007  ! -- Allocate
1008  call mem_allocate(this%nvert, 'NVERT', this%memoryPath)
1009  !
1010  ! -- Initialize
1011  this%nvert = 0
1012  this%ndim = 1
1013  !

◆ connect()

subroutine disv2dmodule::connect ( class(disv2dtype this)
private

Definition at line 591 of file Disv2d.f90.

592  ! -- dummy
593  class(Disv2dType) :: this
594  ! -- local
595  integer(I4B) :: j
596  integer(I4B) :: noder, nrsize
597  integer(I4B) :: narea_eq_zero
598  integer(I4B) :: narea_lt_zero
599  real(DP) :: area
600  !
601  ! -- Initialize
602  narea_eq_zero = 0
603  narea_lt_zero = 0
604  !
605  ! -- Assign the cell area
606  do j = 1, this%nodesuser
607  area = this%get_cell2d_area(j)
608  noder = this%get_nodenumber(j, 0)
609  if (noder > 0) this%area(noder) = area
610  if (area < dzero) then
611  narea_lt_zero = narea_lt_zero + 1
612  write (errmsg, '(a,i0,a)') &
613  &'Calculated CELL2D area less than zero for cell ', j, '.'
614  call store_error(errmsg)
615  end if
616  if (area == dzero) then
617  narea_eq_zero = narea_eq_zero + 1
618  write (errmsg, '(a,i0,a)') &
619  'Calculated CELL2D area is zero for cell ', j, '.'
620  call store_error(errmsg)
621  end if
622  end do
623  !
624  ! -- check for errors
625  if (count_errors() > 0) then
626  if (narea_lt_zero > 0) then
627  write (errmsg, '(i0,a)') narea_lt_zero, &
628  ' cell(s) have an area less than zero. Calculated cell &
629  &areas must be greater than zero. Negative areas often &
630  &mean vertices are not listed in clockwise order.'
631  call store_error(errmsg)
632  end if
633  if (narea_eq_zero > 0) then
634  write (errmsg, '(i0,a)') narea_eq_zero, &
635  ' cell(s) have an area equal to zero. Calculated cell &
636  &areas must be greater than zero. Calculated cell &
637  &areas equal to zero indicate that the cell is not defined &
638  &by a valid polygon.'
639  call store_error(errmsg)
640  end if
641  call store_error_filename(this%input_fname)
642  end if
643  !
644  ! -- create and fill the connections object
645  nrsize = 0
646  if (this%nodes < this%nodesuser) nrsize = this%nodes
647  allocate (this%con)
648  call this%con%disvconnections(this%name_model, this%nodes, &
649  this%nodesuser, 1, nrsize, &
650  this%nvert, this%vertices, this%iavert, &
651  this%javert, this%cellxy, &
652  this%bot, this%bot, &
653  this%nodereduced, this%nodeuser)
654  this%nja = this%con%nja
655  this%njas = this%con%njas
656  !
Here is the call graph for this function:

◆ connection_normal()

subroutine disv2dmodule::connection_normal ( class(disv2dtype this,
integer(i4b), intent(in)  noden,
integer(i4b), intent(in)  nodem,
integer(i4b), intent(in)  ihc,
real(dp), intent(inout)  xcomp,
real(dp), intent(inout)  ycomp,
real(dp), intent(inout)  zcomp,
integer(i4b), intent(in)  ipos 
)
private

The normal points outward from the shared face between noden and nodem.

Parameters
[in]nodencell (reduced nn)
[in]nodemneighbor (reduced nn)
[in]ihchorizontal connection flag

Definition at line 893 of file Disv2d.f90.

895  ! -- dummy
896  class(Disv2dType) :: this
897  integer(I4B), intent(in) :: noden !< cell (reduced nn)
898  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
899  integer(I4B), intent(in) :: ihc !< horizontal connection flag
900  real(DP), intent(inout) :: xcomp
901  real(DP), intent(inout) :: ycomp
902  real(DP), intent(inout) :: zcomp
903  integer(I4B), intent(in) :: ipos
904  ! -- local
905  real(DP) :: angle, dmult
906  !
907  ! -- Set vector components based on ihc
908  if (ihc == 0) then
909  xcomp = dzero
910  ycomp = dzero
911  if (nodem < noden) then
912  !
913  ! -- nodem must be above noden, so upward connection
914  zcomp = done
915  else
916  !
917  ! -- nodem must be below noden, so downward connection
918  zcomp = -done
919  end if
920  else
921  ! -- find from anglex, since anglex is symmetric, need to flip vector
922  ! for lower triangle (nodem < noden)
923  !ipos = this%con%getjaindex(noden, nodem)
924  angle = this%con%anglex(this%con%jas(ipos))
925  dmult = done
926  if (nodem < noden) dmult = -done
927  xcomp = cos(angle) * dmult
928  ycomp = sin(angle) * dmult
929  zcomp = dzero
930  end if
931  !

◆ connection_vector()

subroutine disv2dmodule::connection_vector ( class(disv2dtype this,
integer(i4b), intent(in)  noden,
integer(i4b), intent(in)  nodem,
logical, intent(in)  nozee,
real(dp), intent(in)  satn,
real(dp), intent(in)  satm,
integer(i4b), intent(in)  ihc,
real(dp), intent(inout)  xcomp,
real(dp), intent(inout)  ycomp,
real(dp), intent(inout)  zcomp,
real(dp), intent(inout)  conlen 
)
private

Saturation must be provided to compute cell center vertical coordinates. Also return the straight-line connection length.

Parameters
[in]nodencell (reduced nn)
[in]nodemneighbor (reduced nn)
[in]nozeedo not use z in calculations
[in]satnnot used for disv1d
[in]satmnot used for disv1d
[in]ihchorizontal connection flag
[in,out]xcompx component of connection vector
[in,out]ycompy component of connection vector
[in,out]zcompz component of connection vector
[in,out]conlencalculated straight-line distance between cell centers

Definition at line 939 of file Disv2d.f90.

941  ! -- dummy
942  class(Disv2dType) :: this
943  integer(I4B), intent(in) :: noden !< cell (reduced nn)
944  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
945  logical, intent(in) :: nozee !< do not use z in calculations
946  real(DP), intent(in) :: satn !< not used for disv1d
947  real(DP), intent(in) :: satm !< not used for disv1d
948  integer(I4B), intent(in) :: ihc !< horizontal connection flag
949  real(DP), intent(inout) :: xcomp !< x component of connection vector
950  real(DP), intent(inout) :: ycomp !< y component of connection vector
951  real(DP), intent(inout) :: zcomp !< z component of connection vector
952  real(DP), intent(inout) :: conlen !< calculated straight-line distance between cell centers
953  ! -- local
954  integer(I4B) :: nodeun, nodeum
955  real(DP) :: xn, xm, yn, ym, zn, zm
956 
957  ! horizontal connection, with possible z component due to cell offsets
958  ! and/or water table conditions
959  if (nozee) then
960  zn = dzero
961  zm = dzero
962  else
963  zn = this%bot(noden)
964  zm = this%bot(nodem)
965  end if
966  nodeun = this%get_nodeuser(noden)
967  nodeum = this%get_nodeuser(nodem)
968  xn = this%cellxy(1, nodeun)
969  yn = this%cellxy(2, nodeun)
970  xm = this%cellxy(1, nodeum)
971  ym = this%cellxy(2, nodeum)
972  call line_unit_vector(xn, yn, zn, xm, ym, zm, xcomp, ycomp, zcomp, &
973  conlen)
974 
Here is the call graph for this function:

◆ define_cellverts()

subroutine disv2dmodule::define_cellverts ( class(disv2dtype this,
integer(i4b), dimension(:), intent(in), pointer, contiguous  icell2d,
integer(i4b), dimension(:), intent(in), pointer, contiguous  ncvert,
integer(i4b), dimension(:), intent(in), pointer, contiguous  icvert 
)
private

Definition at line 502 of file Disv2d.f90.

503  ! -- modules
504  use sparsemodule, only: sparsematrix
505  ! -- dummy
506  class(Disv2dType) :: this
507  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: icell2d
508  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: ncvert
509  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: icvert
510  ! -- locals
511  type(sparsematrix) :: vert_spm
512  integer(I4B) :: i, j, ierr
513  integer(I4B) :: icv_idx, startvert, maxnnz = 5
514  !
515  ! -- initialize sparse matrix
516  call vert_spm%init(this%nodes, this%nvert, maxnnz)
517  !
518  ! -- add sparse matrix connections from input memory paths
519  icv_idx = 1
520  do i = 1, this%nodes
521  if (icell2d(i) /= i) call store_error('ICELL2D input sequence violation.')
522  do j = 1, ncvert(i)
523  call vert_spm%addconnection(i, icvert(icv_idx), 0)
524  if (j == 1) then
525  startvert = icvert(icv_idx)
526  elseif (j == ncvert(i) .and. (icvert(icv_idx) /= startvert)) then
527  call vert_spm%addconnection(i, startvert, 0)
528  end if
529  icv_idx = icv_idx + 1
530  end do
531  end do
532  !
533  ! -- allocate and fill iavert and javert
534  call mem_allocate(this%iavert, this%nodes + 1, 'IAVERT', this%memoryPath)
535  call mem_allocate(this%javert, vert_spm%nnz, 'JAVERT', this%memoryPath)
536  call vert_spm%filliaja(this%iavert, this%javert, ierr)
537  call vert_spm%destroy()
538  !
Here is the call graph for this function:

◆ disv2d_cr()

subroutine, public disv2dmodule::disv2d_cr ( class(disbasetype), pointer  dis,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  input_mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Definition at line 94 of file Disv2d.f90.

95  ! -- dummy
96  class(DisBaseType), pointer :: dis
97  character(len=*), intent(in) :: name_model
98  character(len=*), intent(in) :: input_mempath
99  integer(I4B), intent(in) :: inunit
100  integer(I4B), intent(in) :: iout
101  ! -- local
102  type(Disv2dType), pointer :: disnew
103  ! -- formats
104  character(len=*), parameter :: fmtheader = &
105  "(1X, /1X, 'DISV -- VERTEX GRID DISCRETIZATION PACKAGE,', &
106  &' VERSION 1 : 12/23/2015 - INPUT READ FROM MEMPATH: ', A, //)"
107  !
108  allocate (disnew)
109  dis => disnew
110  call disnew%allocate_scalars(name_model, input_mempath)
111  dis%inunit = inunit
112  dis%iout = iout
113  !
114  ! -- If disv enabled
115  if (inunit > 0) then
116  !
117  ! -- Identify package
118  if (iout > 0) then
119  write (iout, fmtheader) dis%input_mempath
120  end if
121  !
122  ! -- load disv
123  call disnew%disv2d_load()
124  end if
125  !
Here is the caller graph for this function:

◆ disv2d_da()

subroutine disv2dmodule::disv2d_da ( class(disv2dtype this)
private

Definition at line 153 of file Disv2d.f90.

154  ! -- modules
158  ! -- dummy
159  class(Disv2dType) :: this
160  ! -- local
161 
162  ! -- Deallocate idm memory
163  call memorystore_remove(this%name_model, 'DISV2D', idm_context)
164 
165  ! -- scalars
166  call mem_deallocate(this%nvert)
167 
168  ! -- arrays
169  call mem_deallocate(this%nodeuser)
170  call mem_deallocate(this%nodereduced)
171  call mem_deallocate(this%bottom)
172  call mem_deallocate(this%idomain)
173 
174  ! -- cdl hack for arrays for vertices and cell2d blocks
175  call mem_deallocate(this%vertices)
176  call mem_deallocate(this%cellxy)
177  call mem_deallocate(this%iavert)
178  call mem_deallocate(this%javert)
179  !
180  ! -- DisBaseType deallocate
181  call this%DisBaseType%dis_da()
subroutine, public memorystore_remove(component, subcomponent, context)
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
Here is the call graph for this function:

◆ disv2d_df()

subroutine disv2dmodule::disv2d_df ( class(disv2dtype this)
private

Definition at line 145 of file Disv2d.f90.

146  ! -- dummy
147  class(Disv2dType) :: this
148  !
149  call this%grid_finalize()
150  !

◆ disv2d_load()

subroutine disv2dmodule::disv2d_load ( class(disv2dtype this)
private

Definition at line 130 of file Disv2d.f90.

131  ! -- dummy
132  class(Disv2dType) :: this
133  !
134  ! -- source input data
135  call this%source_options()
136  call this%source_dimensions()
137  call this%source_griddata()
138  call this%source_vertices()
139  call this%source_cell2d()
140  !

◆ get_cell2d_area()

real(dp) function disv2dmodule::get_cell2d_area ( class(disv2dtype this,
integer(i4b), intent(in)  icell2d 
)
private

A negative result means points are in counter-clockwise orientation. Area is computed from the formula: a = 1/2 *[(x1*y2 + x2*y3 + x3*y4 + ... + xn*y1) - (x2*y1 + x3*y2 + x4*y3 + ... + x1*yn)]

Definition at line 1047 of file Disv2d.f90.

1048  ! -- dummy
1049  class(Disv2dType) :: this
1050  integer(I4B), intent(in) :: icell2d
1051  ! -- return
1052  real(DP) :: area
1053  ! -- local
1054  integer(I4B) :: ivert
1055  integer(I4B) :: nvert
1056  integer(I4B) :: icount
1057  integer(I4B) :: iv1
1058  real(DP) :: x
1059  real(DP) :: y
1060  real(DP) :: x1
1061  real(DP) :: y1
1062  !
1063  area = dzero
1064  nvert = this%iavert(icell2d + 1) - this%iavert(icell2d)
1065  icount = 1
1066  iv1 = this%javert(this%iavert(icell2d))
1067  x1 = this%vertices(1, iv1)
1068  y1 = this%vertices(2, iv1)
1069  do ivert = this%iavert(icell2d), this%iavert(icell2d + 1) - 1
1070  x = this%vertices(1, this%javert(ivert))
1071  if (icount < nvert) then
1072  y = this%vertices(2, this%javert(ivert + 1))
1073  else
1074  y = this%vertices(2, this%javert(this%iavert(icell2d)))
1075  end if
1076  area = area + (x - x1) * (y - y1)
1077  icount = icount + 1
1078  end do
1079  !
1080  icount = 1
1081  do ivert = this%iavert(icell2d), this%iavert(icell2d + 1) - 1
1082  y = this%vertices(2, this%javert(ivert))
1083  if (icount < nvert) then
1084  x = this%vertices(1, this%javert(ivert + 1))
1085  else
1086  x = this%vertices(1, this%javert(this%iavert(icell2d)))
1087  end if
1088  area = area - (x - x1) * (y - y1)
1089  icount = icount + 1
1090  end do
1091  !
1092  area = -done * area * dhalf
1093  !

◆ get_dis_enum()

integer(i4b) function disv2dmodule::get_dis_enum ( class(disv2dtype), intent(in)  this)
private

Definition at line 989 of file Disv2d.f90.

990  use constantsmodule, only: disv2d
991  class(Disv2dType), intent(in) :: this
992  integer(I4B) :: dis_enum
993  dis_enum = disv2d
This module contains simulation constants.
Definition: Constants.f90:9
@ disv2d
DISV2D6 discretization.
Definition: Constants.f90:164

◆ get_dis_type()

subroutine disv2dmodule::get_dis_type ( class(disv2dtype), intent(in)  this,
character(len=*), intent(out)  dis_type 
)
private

Definition at line 979 of file Disv2d.f90.

980  ! -- dummy
981  class(Disv2dType), intent(in) :: this
982  character(len=*), intent(out) :: dis_type
983  !
984  dis_type = "DISV2D"
985  !

◆ get_nodenumber_idx1()

integer(i4b) function disv2dmodule::get_nodenumber_idx1 ( class(disv2dtype), intent(in)  this,
integer(i4b), intent(in)  nodeu,
integer(i4b), intent(in)  icheck 
)
private

Definition at line 859 of file Disv2d.f90.

860  ! return
861  integer(I4B) :: nodenumber
862  ! dummy
863  class(Disv2dType), intent(in) :: this
864  integer(I4B), intent(in) :: nodeu
865  integer(I4B), intent(in) :: icheck
866  ! local
867 
868  ! check the node number if requested
869  if (icheck /= 0) then
870 
871  ! If within valid range, convert to reduced nodenumber
872  if (nodeu < 1 .or. nodeu > this%nodesuser) then
873  nodenumber = 0
874  write (errmsg, '(a,i0,a,i0,a)') &
875  'Node number (', nodeu, ') is less than 1 or greater than nodes (', &
876  this%nodesuser, ').'
877  call store_error(errmsg)
878  else
879  nodenumber = nodeu
880  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
881  end if
882  else
883  nodenumber = nodeu
884  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
885  end if
886 
Here is the call graph for this function:

◆ get_polyverts()

subroutine disv2dmodule::get_polyverts ( class(disv2dtype), intent(inout)  this,
integer(i4b), intent(in)  ic,
real(dp), dimension(:, :), intent(out), allocatable  polyverts,
logical(lgp), intent(in), optional  closed 
)
private
Parameters
[in]iccell number (reduced)
[out]polyvertspolygon vertices (column-major indexing)
[in]closedwhether to close the polygon, duplicating a vertex (default false)

Definition at line 1253 of file Disv2d.f90.

1254  ! -- dummy
1255  class(Disv2dType), intent(inout) :: this
1256  integer(I4B), intent(in) :: ic !< cell number (reduced)
1257  real(DP), allocatable, intent(out) :: polyverts(:, :) !< polygon vertices (column-major indexing)
1258  logical(LGP), intent(in), optional :: closed !< whether to close the polygon, duplicating a vertex (default false)
1259  ! -- local
1260  integer(I4B) :: icu, icu2d, iavert, nverts, m, j
1261  logical(LGP) :: lclosed
1262  !
1263  ! count vertices
1264  icu = this%get_nodeuser(ic)
1265  icu2d = icu - ((icu - 1) / this%nodes) * this%nodes
1266  nverts = this%iavert(icu2d + 1) - this%iavert(icu2d) - 1
1267  if (nverts .le. 0) nverts = nverts + size(this%javert)
1268  !
1269  ! check closed option
1270  if (.not. (present(closed))) then
1271  lclosed = .false.
1272  else
1273  lclosed = closed
1274  end if
1275  !
1276  ! allocate vertices array
1277  if (lclosed) then
1278  allocate (polyverts(2, nverts + 1))
1279  else
1280  allocate (polyverts(2, nverts))
1281  end if
1282  !
1283  ! set vertices
1284  iavert = this%iavert(icu2d)
1285  do m = 1, nverts
1286  j = this%javert(iavert - 1 + m)
1287  polyverts(:, m) = (/this%vertices(1, j), this%vertices(2, j)/)
1288  end do
1289  !
1290  ! close if enabled
1291  if (lclosed) &
1292  polyverts(:, nverts + 1) = polyverts(:, 1)
1293  !

◆ grid_finalize()

subroutine disv2dmodule::grid_finalize ( class(disv2dtype this)
private

Definition at line 386 of file Disv2d.f90.

387  ! dummy
388  class(Disv2dType) :: this
389  ! locals
390  integer(I4B) :: node, noder, j, ncell_count
391  ! formats
392  character(len=*), parameter :: fmtnr = &
393  "(/1x, 'The specified IDOMAIN results in a reduced number of cells.',&
394  &/1x, 'Number of user nodes: ',I0,&
395  &/1X, 'Number of nodes in solution: ', I0, //)"
396 
397  ! count active cells and set nodes to that number
398  ncell_count = 0
399  do j = 1, this%nodesuser
400  if (this%idomain(j) > 0) ncell_count = ncell_count + 1
401  end do
402  this%nodes = ncell_count
403 
404  ! Check to make sure nodes is a valid number
405  if (ncell_count == 0) then
406  call store_error('Model does not have any active nodes. &
407  &Ensure IDOMAIN array has some values greater &
408  &than zero.')
409  call store_error_filename(this%input_fname)
410  end if
411 
412  ! Write message if reduced grid
413  if (this%nodes < this%nodesuser) then
414  write (this%iout, fmtnr) this%nodesuser, this%nodes
415  end if
416 
417  ! Array size is now known, so allocate
418  call this%allocate_arrays()
419 
420  ! Fill the nodereduced array with the reduced nodenumber, or
421  ! a negative number to indicate it is a pass-through cell, or
422  ! a zero to indicate that the cell is excluded from the
423  ! solution.
424  if (this%nodes < this%nodesuser) then
425  node = 1
426  noder = 1
427  do j = 1, this%nodesuser
428  if (this%idomain(j) > 0) then
429  this%nodereduced(node) = noder
430  noder = noder + 1
431  else
432  this%nodereduced(node) = 0
433  end if
434  node = node + 1
435  end do
436  end if
437 
438  ! allocate and fill nodeuser if a reduced grid
439  if (this%nodes < this%nodesuser) then
440  node = 1
441  noder = 1
442  do j = 1, this%nodesuser
443  if (this%idomain(j) > 0) then
444  this%nodeuser(noder) = node
445  noder = noder + 1
446  end if
447  node = node + 1
448  end do
449  end if
450 
451  ! Move bottom into bot
452  ! and set x and y center coordinates
453  node = 0
454  do j = 1, this%nodesuser
455  node = node + 1
456  noder = node
457  if (this%nodes < this%nodesuser) noder = this%nodereduced(node)
458  if (noder <= 0) cycle
459  this%bot(noder) = this%bottom(j)
460  this%xc(noder) = this%cellxy(1, j)
461  this%yc(noder) = this%cellxy(2, j)
462  end do
463 
464  ! Build connections
465  call this%connect()
466 
Here is the call graph for this function:

◆ log_dimensions()

subroutine disv2dmodule::log_dimensions ( class(disv2dtype this,
type(disvfoundtype), intent(in)  found 
)
private

Definition at line 325 of file Disv2d.f90.

326  ! -- dummy
327  class(Disv2dType) :: this
328  type(DisvFoundType), intent(in) :: found
329  !
330  write (this%iout, '(1x,a)') 'Setting Discretization Dimensions'
331  !
332  if (found%nodes) then
333  write (this%iout, '(4x,a,i0)') 'NODES = ', this%nodesuser
334  end if
335  !
336  if (found%nvert) then
337  write (this%iout, '(4x,a,i0)') 'NVERT = ', this%nvert
338  end if
339  !
340  write (this%iout, '(1x,a,/)') 'End Setting Discretization Dimensions'
341  !

◆ log_griddata()

subroutine disv2dmodule::log_griddata ( class(disv2dtype this,
type(disvfoundtype), intent(in)  found 
)
private

Definition at line 365 of file Disv2d.f90.

366  ! -- dummy
367  class(Disv2dType) :: this
368  type(DisvFoundType), intent(in) :: found
369  !
370  write (this%iout, '(1x,a)') 'Setting Discretization Griddata'
371  !
372  if (found%bottom) then
373  write (this%iout, '(4x,a)') 'BOTTOM set from input file'
374  end if
375  !
376  if (found%idomain) then
377  write (this%iout, '(4x,a)') 'IDOMAIN set from input file'
378  end if
379  !
380  write (this%iout, '(1x,a,/)') 'End Setting Discretization Griddata'
381  !

◆ log_options()

subroutine disv2dmodule::log_options ( class(disv2dtype this,
type(disvfoundtype), intent(in)  found 
)
private

Definition at line 240 of file Disv2d.f90.

241  ! -- dummy
242  class(Disv2dType) :: this
243  type(DisvFoundType), intent(in) :: found
244  !
245  write (this%iout, '(1x,a)') 'Setting Discretization Options'
246  !
247  if (found%length_units) then
248  write (this%iout, '(4x,a,i0)') 'Model length unit [0=UND, 1=FEET, &
249  &2=METERS, 3=CENTIMETERS] set as ', this%lenuni
250  end if
251  !
252  if (found%nogrb) then
253  write (this%iout, '(4x,a,i0)') 'Binary grid file [0=GRB, 1=NOGRB] &
254  &set as ', this%nogrb
255  end if
256  !
257  if (found%xorigin) then
258  write (this%iout, '(4x,a,G0)') 'XORIGIN = ', this%xorigin
259  end if
260  !
261  if (found%yorigin) then
262  write (this%iout, '(4x,a,G0)') 'YORIGIN = ', this%yorigin
263  end if
264  !
265  if (found%angrot) then
266  write (this%iout, '(4x,a,G0)') 'ANGROT = ', this%angrot
267  end if
268  !
269  write (this%iout, '(1x,a,/)') 'End Setting Discretization Options'
270  !

◆ nodeu_from_cellid()

integer(i4b) function disv2dmodule::nodeu_from_cellid ( class(disv2dtype this,
character(len=*), intent(inout)  cellid,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
logical, intent(in), optional  flag_string,
logical, intent(in), optional  allow_zero 
)
private

If flag_string is present and true, the first token may be non-numeric (e.g. boundary name). In this case, return -2.

If allow_zero is present and true, and all indices are zero, the result can be zero. If allow_zero is false, a zero in any index is an error.

Definition at line 1180 of file Disv2d.f90.

1182  ! -- return
1183  integer(I4B) :: nodeu
1184  ! -- dummy
1185  class(Disv2dType) :: this
1186  character(len=*), intent(inout) :: cellid
1187  integer(I4B), intent(in) :: inunit
1188  integer(I4B), intent(in) :: iout
1189  logical, optional, intent(in) :: flag_string
1190  logical, optional, intent(in) :: allow_zero
1191  ! -- local
1192  integer(I4B) :: j, nodes
1193  integer(I4B) :: lloclocal, ndum, istat, n
1194  integer(I4B) :: istart, istop
1195  real(DP) :: r
1196  !
1197  if (present(flag_string)) then
1198  if (flag_string) then
1199  ! Check to see if first token in cellid can be read as an integer.
1200  lloclocal = 1
1201  call urword(cellid, lloclocal, istart, istop, 1, ndum, r, iout, inunit)
1202  read (cellid(istart:istop), *, iostat=istat) n
1203  if (istat /= 0) then
1204  ! First token in cellid is not an integer; return flag to this effect.
1205  nodeu = -2
1206  return
1207  end if
1208  end if
1209  end if
1210  !
1211  nodes = this%mshape(1)
1212  !
1213  lloclocal = 1
1214  call urword(cellid, lloclocal, istart, istop, 2, j, r, iout, inunit)
1215  !
1216  if (j == 0) then
1217  if (present(allow_zero)) then
1218  if (allow_zero) then
1219  nodeu = 0
1220  return
1221  end if
1222  end if
1223  end if
1224  !
1225  errmsg = ''
1226  !
1227  if (j < 1 .or. j > nodes) then
1228  write (errmsg, '(a,1x,a,i0,a)') &
1229  trim(adjustl(errmsg)), 'Cell2d number in list (', j, &
1230  ') is outside of the grid.'
1231  end if
1232  !
1233  nodeu = get_node(1, 1, j, 1, 1, nodes)
1234  !
1235  if (nodeu < 1 .or. nodeu > this%nodesuser) then
1236  write (errmsg, '(a,1x,a,i0,a)') &
1237  trim(adjustl(errmsg)), &
1238  "Cell number cannot be determined for cellid ("// &
1239  trim(adjustl(cellid))//") and results in a user "// &
1240  "node number (", nodeu, ") that is outside of the grid."
1241  end if
1242  !
1243  if (len_trim(adjustl(errmsg)) > 0) then
1244  call store_error(errmsg)
1245  call store_error_unit(inunit)
1246  end if
1247  !
Here is the call graph for this function:

◆ nodeu_from_string()

integer(i4b) function disv2dmodule::nodeu_from_string ( class(disv2dtype this,
integer(i4b), intent(inout)  lloc,
integer(i4b), intent(inout)  istart,
integer(i4b), intent(inout)  istop,
integer(i4b), intent(in)  in,
integer(i4b), intent(in)  iout,
character(len=*), intent(inout)  line,
logical, intent(in), optional  flag_string,
logical, intent(in), optional  allow_zero 
)
private

Parse layer and within-layer cell number and return user nodenumber. If flag_string is present and true, the first token may be non-numeric (e.g. boundary name). In this case, return -2.

Definition at line 1102 of file Disv2d.f90.

1104  ! -- dummy
1105  class(Disv2dType) :: this
1106  integer(I4B), intent(inout) :: lloc
1107  integer(I4B), intent(inout) :: istart
1108  integer(I4B), intent(inout) :: istop
1109  integer(I4B), intent(in) :: in
1110  integer(I4B), intent(in) :: iout
1111  character(len=*), intent(inout) :: line
1112  logical, optional, intent(in) :: flag_string
1113  logical, optional, intent(in) :: allow_zero
1114  integer(I4B) :: nodeu
1115  ! -- local
1116  integer(I4B) :: j, nodes
1117  integer(I4B) :: lloclocal, ndum, istat, n
1118  real(DP) :: r
1119  !
1120  if (present(flag_string)) then
1121  if (flag_string) then
1122  ! Check to see if first token in line can be read as an integer.
1123  lloclocal = lloc
1124  call urword(line, lloclocal, istart, istop, 1, ndum, r, iout, in)
1125  read (line(istart:istop), *, iostat=istat) n
1126  if (istat /= 0) then
1127  ! First token in line is not an integer; return flag to this effect.
1128  nodeu = -2
1129  return
1130  end if
1131  end if
1132  end if
1133  !
1134  nodes = this%mshape(1)
1135  !
1136  call urword(line, lloc, istart, istop, 2, j, r, iout, in)
1137  !
1138  if (j == 0) then
1139  if (present(allow_zero)) then
1140  if (allow_zero) then
1141  nodeu = 0
1142  return
1143  end if
1144  end if
1145  end if
1146  !
1147  errmsg = ''
1148  !
1149  if (j < 1 .or. j > nodes) then
1150  write (errmsg, '(a,1x,a,i0,a)') &
1151  trim(adjustl(errmsg)), 'Cell number in list (', j, &
1152  ') is outside of the grid.'
1153  end if
1154  !
1155  nodeu = get_node(1, 1, j, 1, 1, nodes)
1156  !
1157  if (nodeu < 1 .or. nodeu > this%nodesuser) then
1158  write (errmsg, '(a,1x,a,i0,a)') &
1159  trim(adjustl(errmsg)), &
1160  "Node number in list (", nodeu, ") is outside of the grid. "// &
1161  "Cell number cannot be determined in line '"// &
1162  trim(adjustl(line))//"'."
1163  end if
1164  !
1165  if (len_trim(adjustl(errmsg)) > 0) then
1166  call store_error(errmsg)
1167  call store_error_unit(in)
1168  end if
1169  !
Here is the call graph for this function:

◆ nodeu_to_array()

subroutine disv2dmodule::nodeu_to_array ( class(disv2dtype this,
integer(i4b), intent(in)  nodeu,
integer(i4b), dimension(:), intent(inout)  arr 
)
private

Definition at line 831 of file Disv2d.f90.

832  ! -- dummy
833  class(Disv2dType) :: this
834  integer(I4B), intent(in) :: nodeu
835  integer(I4B), dimension(:), intent(inout) :: arr
836  ! -- local
837  integer(I4B) :: isize
838  integer(I4B) :: i, j, k
839  !
840  ! -- check the size of arr
841  isize = size(arr)
842  if (isize /= this%ndim) then
843  write (errmsg, '(a,i0,a,i0,a)') &
844  'Program error: nodeu_to_array size of array (', isize, &
845  ') is not equal to the discretization dimension (', this%ndim, ').'
846  call store_error(errmsg, terminate=.true.)
847  end if
848  !
849  ! -- get k, i, j
850  call get_ijk(nodeu, 1, this%nodes, 1, i, j, k)
851  !
852  ! -- fill array
853  arr(1) = j
854  !
Here is the call graph for this function:

◆ nodeu_to_string()

subroutine disv2dmodule::nodeu_to_string ( class(disv2dtype this,
integer(i4b), intent(in)  nodeu,
character(len=*), intent(inout)  str 
)

Definition at line 814 of file Disv2d.f90.

815  ! -- dummy
816  class(Disv2dType) :: this
817  integer(I4B), intent(in) :: nodeu
818  character(len=*), intent(inout) :: str
819  ! -- local
820  integer(I4B) :: i, j, k
821  character(len=10) :: jstr
822  !
823  call get_ijk(nodeu, 1, this%nodes, 1, i, j, k)
824  write (jstr, '(i10)') j
825  str = '('//trim(adjustl(jstr))//')'
826  !
Here is the call graph for this function:

◆ record_array()

subroutine disv2dmodule::record_array ( class(disv2dtype), intent(inout)  this,
real(dp), dimension(:), intent(inout), pointer, contiguous  darray,
integer(i4b), intent(in)  iout,
integer(i4b), intent(in)  iprint,
integer(i4b), intent(in)  idataun,
character(len=*), intent(in)  aname,
character(len=*), intent(in)  cdatafmp,
integer(i4b), intent(in)  nvaluesp,
integer(i4b), intent(in)  nwidthp,
character(len=*), intent(in)  editdesc,
real(dp), intent(in)  dinact 
)
private

The array is written to a formatted or unformatted external file depending on the arguments.

Parameters
[in,out]darraydouble precision array to record
[in]ioutascii output unit number
[in]iprintwhether to print the array
[in]idataunbinary output unit number, if negative don't write by layers, write entire array
[in]anametext descriptor
[in]cdatafmpwrite format
[in]nvaluespvalues per line
[in]nwidthpnumber width
[in]editdescformat type (I, G, F, S, E)
[in]dinactdouble precision value for cells excluded from model domain

Definition at line 1301 of file Disv2d.f90.

1303  ! -- dummy
1304  class(Disv2dType), intent(inout) :: this
1305  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray !< double precision array to record
1306  integer(I4B), intent(in) :: iout !< ascii output unit number
1307  integer(I4B), intent(in) :: iprint !< whether to print the array
1308  integer(I4B), intent(in) :: idataun !< binary output unit number, if negative don't write by layers, write entire array
1309  character(len=*), intent(in) :: aname !< text descriptor
1310  character(len=*), intent(in) :: cdatafmp !< write format
1311  integer(I4B), intent(in) :: nvaluesp !< values per line
1312  integer(I4B), intent(in) :: nwidthp !< number width
1313  character(len=*), intent(in) :: editdesc !< format type (I, G, F, S, E)
1314  real(DP), intent(in) :: dinact !< double precision value for cells excluded from model domain
1315  ! -- local
1316  integer(I4B) :: k, ifirst
1317  integer(I4B) :: nlay
1318  integer(I4B) :: nrow
1319  integer(I4B) :: ncol
1320  integer(I4B) :: nval
1321  integer(I4B) :: nodeu, noder
1322  integer(I4B) :: istart, istop
1323  real(DP), dimension(:), pointer, contiguous :: dtemp
1324  ! -- formats
1325  character(len=*), parameter :: fmthsv = &
1326  "(1X,/1X,a,' WILL BE SAVED ON UNIT ',I4, &
1327  &' AT END OF TIME STEP',I5,', STRESS PERIOD ',I4)"
1328  !
1329  ! -- set variables
1330  nlay = 1
1331  nrow = 1
1332  ncol = this%mshape(1)
1333  !
1334  ! -- If this is a reduced model, then copy the values from darray into
1335  ! dtemp.
1336  if (this%nodes < this%nodesuser) then
1337  nval = this%nodes
1338  dtemp => this%dbuff
1339  do nodeu = 1, this%nodesuser
1340  noder = this%get_nodenumber(nodeu, 0)
1341  if (noder <= 0) then
1342  dtemp(nodeu) = dinact
1343  cycle
1344  end if
1345  dtemp(nodeu) = darray(noder)
1346  end do
1347  else
1348  nval = this%nodes
1349  dtemp => darray
1350  end if
1351  !
1352  ! -- Print to iout if iprint /= 0
1353  if (iprint /= 0) then
1354  istart = 1
1355  do k = 1, nlay
1356  istop = istart + nrow * ncol - 1
1357  call ulaprufw(ncol, nrow, kstp, kper, k, iout, dtemp(istart:istop), &
1358  aname, cdatafmp, nvaluesp, nwidthp, editdesc)
1359  istart = istop + 1
1360  end do
1361  end if
1362  !
1363  ! -- Save array to an external file.
1364  if (idataun > 0) then
1365  ! -- write to binary file by layer
1366  ifirst = 1
1367  istart = 1
1368  do k = 1, nlay
1369  istop = istart + nrow * ncol - 1
1370  if (ifirst == 1) write (iout, fmthsv) &
1371  trim(adjustl(aname)), idataun, &
1372  kstp, kper
1373  ifirst = 0
1374  call ulasav(dtemp(istart:istop), aname, kstp, kper, &
1375  pertim, totim, ncol, nrow, k, idataun)
1376  istart = istop + 1
1377  end do
1378  elseif (idataun < 0) then
1379  !
1380  ! -- write entire array as one record
1381  call ubdsv1(kstp, kper, aname, -idataun, dtemp, ncol, nrow, nlay, &
1382  iout, delt, pertim, totim)
1383  end if
1384  !
Here is the call graph for this function:

◆ record_srcdst_list_header()

subroutine disv2dmodule::record_srcdst_list_header ( class(disv2dtype this,
character(len=16), intent(in)  text,
character(len=16), intent(in)  textmodel,
character(len=16), intent(in)  textpackage,
character(len=16), intent(in)  dstmodel,
character(len=16), intent(in)  dstpackage,
integer(i4b), intent(in)  naux,
character(len=16), dimension(:), intent(in)  auxtxt,
integer(i4b), intent(in)  ibdchn,
integer(i4b), intent(in)  nlist,
integer(i4b), intent(in)  iout 
)
private

Definition at line 1389 of file Disv2d.f90.

1392  ! -- dummy
1393  class(Disv2dType) :: this
1394  character(len=16), intent(in) :: text
1395  character(len=16), intent(in) :: textmodel
1396  character(len=16), intent(in) :: textpackage
1397  character(len=16), intent(in) :: dstmodel
1398  character(len=16), intent(in) :: dstpackage
1399  integer(I4B), intent(in) :: naux
1400  character(len=16), dimension(:), intent(in) :: auxtxt
1401  integer(I4B), intent(in) :: ibdchn
1402  integer(I4B), intent(in) :: nlist
1403  integer(I4B), intent(in) :: iout
1404  ! -- local
1405  integer(I4B) :: nlay, nrow, ncol
1406  !
1407  nlay = 1
1408  nrow = 1
1409  ncol = this%mshape(1)
1410  !
1411  ! -- Use ubdsv06 to write list header
1412  call ubdsv06(kstp, kper, text, textmodel, textpackage, dstmodel, dstpackage, &
1413  ibdchn, naux, auxtxt, ncol, nrow, nlay, &
1414  nlist, iout, delt, pertim, totim)
1415  !
Here is the call graph for this function:

◆ source_cell2d()

subroutine disv2dmodule::source_cell2d ( class(disv2dtype this)

Definition at line 543 of file Disv2d.f90.

544  ! -- dummy
545  class(Disv2dType) :: this
546  ! -- locals
547  integer(I4B), dimension(:), contiguous, pointer :: icell2d => null()
548  integer(I4B), dimension(:), contiguous, pointer :: ncvert => null()
549  integer(I4B), dimension(:), contiguous, pointer :: icvert => null()
550  real(DP), dimension(:), contiguous, pointer :: cell_x => null()
551  real(DP), dimension(:), contiguous, pointer :: cell_y => null()
552  integer(I4B) :: i
553  !
554  ! -- set pointers to input path ncvert and icvert
555  call mem_setptr(icell2d, 'ICELL2D', this%input_mempath)
556  call mem_setptr(ncvert, 'NCVERT', this%input_mempath)
557  call mem_setptr(icvert, 'ICVERT', this%input_mempath)
558  !
559  ! --
560  if (associated(icell2d) .and. associated(ncvert) &
561  .and. associated(icvert)) then
562  call this%define_cellverts(icell2d, ncvert, icvert)
563  else
564  call store_error('Required cell vertex array(s) [ICELL2D, NCVERT, ICVERT] &
565  &not found.')
566  end if
567  !
568  ! -- copy cell center idm sourced values to local arrays
569  call mem_setptr(cell_x, 'XC', this%input_mempath)
570  call mem_setptr(cell_y, 'YC', this%input_mempath)
571  !
572  ! -- set cell centers
573  if (associated(cell_x) .and. associated(cell_y)) then
574  do i = 1, this%nodesuser
575  this%cellxy(1, i) = cell_x(i)
576  this%cellxy(2, i) = cell_y(i)
577  end do
578  else
579  call store_error('Required cell center arrays not found.')
580  end if
581  !
582  ! -- log
583  if (this%iout > 0) then
584  write (this%iout, '(1x,a)') 'Discretization Cell2d data loaded'
585  end if
586  !
Here is the call graph for this function:

◆ source_dimensions()

subroutine disv2dmodule::source_dimensions ( class(disv2dtype this)
private

Definition at line 275 of file Disv2d.f90.

276  ! -- dummy
277  class(Disv2dType) :: this
278  ! -- locals
279  integer(I4B) :: j
280  type(DisvFoundType) :: found
281  !
282  ! -- update defaults with idm sourced values
283  call mem_set_value(this%nodes, 'NODES', this%input_mempath, found%nodes)
284  call mem_set_value(this%nvert, 'NVERT', this%input_mempath, found%nvert)
285  !
286  ! -- log simulation values
287  if (this%iout > 0) then
288  call this%log_dimensions(found)
289  end if
290  !
291  ! -- verify dimensions were set
292  if (this%nodes < 1) then
293  call store_error( &
294  'NODES was not specified or was specified incorrectly.')
295  call store_error_filename(this%input_fname)
296  end if
297  if (this%nvert < 1) then
298  call store_error( &
299  'NVERT was not specified or was specified incorrectly.')
300  call store_error_filename(this%input_fname)
301  end if
302  !
303  ! -- Calculate nodesuser
304  this%nodesuser = this%nodes
305  !
306  ! -- Allocate non-reduced vectors for disv
307  call mem_allocate(this%idomain, this%nodes, 'IDOMAIN', &
308  this%memoryPath)
309  call mem_allocate(this%bottom, this%nodes, 'BOTTOM', &
310  this%memoryPath)
311  !
312  ! -- Allocate vertices array
313  call mem_allocate(this%vertices, 2, this%nvert, 'VERTICES', this%memoryPath)
314  call mem_allocate(this%cellxy, 2, this%nodesuser, 'CELLXY', this%memoryPath)
315  !
316  ! -- initialize all cells to be active (idomain = 1)
317  do j = 1, this%nodesuser
318  this%idomain(j) = 1
319  end do
320  !
Here is the call graph for this function:

◆ source_griddata()

subroutine disv2dmodule::source_griddata ( class(disv2dtype this)
private

Definition at line 346 of file Disv2d.f90.

347  ! -- dummy
348  class(Disv2dType) :: this
349  ! -- locals
350  type(DisvFoundType) :: found
351  !
352  ! -- update defaults with idm sourced values
353  call mem_set_value(this%bottom, 'BOTTOM', this%input_mempath, found%bottom)
354  call mem_set_value(this%idomain, 'IDOMAIN', this%input_mempath, found%idomain)
355  !
356  ! -- log simulation values
357  if (this%iout > 0) then
358  call this%log_griddata(found)
359  end if
360  !

◆ source_options()

subroutine disv2dmodule::source_options ( class(disv2dtype this)

Definition at line 215 of file Disv2d.f90.

216  ! -- dummy
217  class(Disv2dType) :: this
218  ! -- locals
219  character(len=LENVARNAME), dimension(3) :: lenunits = &
220  &[character(len=LENVARNAME) :: 'FEET', 'METERS', 'CENTIMETERS']
221  type(disvfoundtype) :: found
222  !
223  ! -- update defaults with idm sourced values
224  call mem_set_value(this%lenuni, 'LENGTH_UNITS', this%input_mempath, &
225  lenunits, found%length_units)
226  call mem_set_value(this%nogrb, 'NOGRB', this%input_mempath, found%nogrb)
227  call mem_set_value(this%xorigin, 'XORIGIN', this%input_mempath, found%xorigin)
228  call mem_set_value(this%yorigin, 'YORIGIN', this%input_mempath, found%yorigin)
229  call mem_set_value(this%angrot, 'ANGROT', this%input_mempath, found%angrot)
230  !
231  ! -- log values to list file
232  if (this%iout > 0) then
233  call this%log_options(found)
234  end if
235  !

◆ source_vertices()

subroutine disv2dmodule::source_vertices ( class(disv2dtype this)
private

Definition at line 471 of file Disv2d.f90.

472  ! -- dummy
473  class(Disv2dType) :: this
474  ! -- local
475  integer(I4B) :: i
476  real(DP), dimension(:), contiguous, pointer :: vert_x => null()
477  real(DP), dimension(:), contiguous, pointer :: vert_y => null()
478  !
479  ! -- set pointers to memory manager input arrays
480  call mem_setptr(vert_x, 'XV', this%input_mempath)
481  call mem_setptr(vert_y, 'YV', this%input_mempath)
482  !
483  ! -- set vertices 2d array
484  if (associated(vert_x) .and. associated(vert_y)) then
485  do i = 1, this%nvert
486  this%vertices(1, i) = vert_x(i)
487  this%vertices(2, i) = vert_y(i)
488  end do
489  else
490  call store_error('Required Vertex arrays not found.')
491  end if
492  !
493  ! -- log
494  if (this%iout > 0) then
495  write (this%iout, '(1x,a)') 'Discretization Vertex data loaded'
496  end if
497  !
Here is the call graph for this function:

◆ write_grb()

subroutine disv2dmodule::write_grb ( class(disv2dtype this,
integer(i4b), dimension(:), intent(in)  icelltype 
)
private

Definition at line 661 of file Disv2d.f90.

662  ! -- modules
663  use openspecmodule, only: access, form
664  use constantsmodule, only: lenbigline
665  ! -- dummy
666  class(Disv2dType) :: this
667  integer(I4B), dimension(:), intent(in) :: icelltype
668  ! -- local
669  integer(I4B) :: iunit, i, ntxt, version
670  integer(I4B), parameter :: lentxt = 100
671  character(len=50) :: txthdr
672  character(len=lentxt) :: txt
673  character(len=LINELENGTH) :: fname
674  character(len=LENBIGLINE) :: crs
675  logical(LGP) :: found_crs
676  ! -- formats
677  character(len=*), parameter :: fmtgrdsave = &
678  "(4X,'BINARY GRID INFORMATION WILL BE WRITTEN TO:', &
679  &/,6X,'UNIT NUMBER: ', I0,/,6X, 'FILE NAME: ', A)"
680  !
681  ! -- Initialize
682  version = 1
683  ntxt = 18
684  !
685  call mem_set_value(crs, 'CRS', this%input_mempath, found_crs)
686  !
687  ! -- set version
688  if (found_crs) then
689  ntxt = ntxt + 1
690  version = 2
691  end if
692  !
693  ! -- Open the file
694  fname = trim(this%output_fname)
695  iunit = getunit()
696  write (this%iout, fmtgrdsave) iunit, trim(adjustl(fname))
697  call openfile(iunit, this%iout, trim(adjustl(fname)), 'DATA(BINARY)', &
698  form, access, 'REPLACE')
699  !
700  ! -- write header information
701  write (txthdr, '(a)') 'GRID DISV2D'
702  txthdr(50:50) = new_line('a')
703  write (iunit) txthdr
704  write (txthdr, '(a)') 'VERSION 1'
705  txthdr(50:50) = new_line('a')
706  write (iunit) txthdr
707  write (txthdr, '(a, i0)') 'NTXT ', ntxt
708  txthdr(50:50) = new_line('a')
709  write (iunit) txthdr
710  write (txthdr, '(a, i0)') 'LENTXT ', lentxt
711  txthdr(50:50) = new_line('a')
712  write (iunit) txthdr
713  !
714  ! -- write variable definitions
715  write (txt, '(3a, i0)') 'NCELLS ', 'INTEGER ', 'NDIM 0 # ', this%nodesuser
716  txt(lentxt:lentxt) = new_line('a')
717  write (iunit) txt
718  write (txt, '(3a, i0)') 'NODES ', 'INTEGER ', 'NDIM 0 # ', this%nodes
719  txt(lentxt:lentxt) = new_line('a')
720  write (iunit) txt
721  write (txt, '(3a, i0)') 'NVERT ', 'INTEGER ', 'NDIM 0 # ', this%nvert
722  txt(lentxt:lentxt) = new_line('a')
723  write (iunit) txt
724  write (txt, '(3a, i0)') 'NJAVERT ', 'INTEGER ', 'NDIM 0 # ', size(this%javert)
725  txt(lentxt:lentxt) = new_line('a')
726  write (iunit) txt
727  write (txt, '(3a, i0)') 'NJA ', 'INTEGER ', 'NDIM 0 # ', this%con%nja
728  txt(lentxt:lentxt) = new_line('a')
729  write (iunit) txt
730  write (txt, '(3a, 1pg25.15e3)') &
731  'XORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%xorigin
732  txt(lentxt:lentxt) = new_line('a')
733  write (iunit) txt
734  write (txt, '(3a, 1pg25.15e3)') &
735  'YORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%yorigin
736  txt(lentxt:lentxt) = new_line('a')
737  write (iunit) txt
738  write (txt, '(3a, 1pg25.15e3)') 'ANGROT ', 'DOUBLE ', 'NDIM 0 # ', this%angrot
739  txt(lentxt:lentxt) = new_line('a')
740  write (iunit) txt
741  write (txt, '(3a, i0)') 'BOTM ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
742  txt(lentxt:lentxt) = new_line('a')
743  write (iunit) txt
744  write (txt, '(3a, i0)') 'VERTICES ', 'DOUBLE ', 'NDIM 2 2 ', this%nvert
745  txt(lentxt:lentxt) = new_line('a')
746  write (iunit) txt
747  write (txt, '(3a, i0)') 'CELLX ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
748  txt(lentxt:lentxt) = new_line('a')
749  write (iunit) txt
750  write (txt, '(3a, i0)') 'CELLY ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
751  txt(lentxt:lentxt) = new_line('a')
752  write (iunit) txt
753  write (txt, '(3a, i0)') 'IAVERT ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
754  txt(lentxt:lentxt) = new_line('a')
755  write (iunit) txt
756  write (txt, '(3a, i0)') 'JAVERT ', 'INTEGER ', 'NDIM 1 ', size(this%javert)
757  txt(lentxt:lentxt) = new_line('a')
758  write (iunit) txt
759  write (txt, '(3a, i0)') 'IA ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
760  txt(lentxt:lentxt) = new_line('a')
761  write (iunit) txt
762  write (txt, '(3a, i0)') 'JA ', 'INTEGER ', 'NDIM 1 ', size(this%con%jausr)
763  txt(lentxt:lentxt) = new_line('a')
764  write (iunit) txt
765  write (txt, '(3a, i0)') 'IDOMAIN ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
766  txt(lentxt:lentxt) = new_line('a')
767  write (iunit) txt
768  write (txt, '(3a, i0)') 'ICELLTYPE ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
769  txt(lentxt:lentxt) = new_line('a')
770  write (iunit) txt
771  !
772  ! -- if version 2 write character array headers
773  if (version == 2) then
774  if (found_crs) then
775  write (txt, '(3a, i0)') 'CRS ', 'CHARACTER ', 'NDIM 1 ', &
776  len_trim(crs)
777  txt(lentxt:lentxt) = new_line('a')
778  write (iunit) txt
779  end if
780  end if
781  !
782  ! -- write data
783  write (iunit) this%nodesuser ! ncells
784  write (iunit) this%nodes ! nodes
785  write (iunit) this%nvert ! nvert
786  write (iunit) size(this%javert) ! njavert
787  write (iunit) this%nja ! nja
788  write (iunit) this%xorigin ! xorigin
789  write (iunit) this%yorigin ! yorigin
790  write (iunit) this%angrot ! angrot
791  write (iunit) this%bottom ! botm
792  write (iunit) this%vertices ! vertices
793  write (iunit) (this%cellxy(1, i), i=1, this%nodesuser) ! cellx
794  write (iunit) (this%cellxy(2, i), i=1, this%nodesuser) ! celly
795  write (iunit) this%iavert ! iavert
796  write (iunit) this%javert ! javert
797  write (iunit) this%con%iausr ! iausr
798  write (iunit) this%con%jausr ! jausr
799  write (iunit) this%idomain ! idomain
800  write (iunit) icelltype ! icelltype
801  !
802  ! -- if version 2 write character array data
803  if (version == 2) then
804  if (found_crs) write (iunit) trim(crs) ! crs user input
805  end if
806  !
807  ! -- Close the file
808  close (iunit)
809  !
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Here is the call graph for this function: