MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
CellWithNbrs.f90
Go to the documentation of this file.
2  use kindmodule, only: i4b, lgp
5  implicit none
6  private
7 
8  integer(I4B), parameter :: defaultcapacity = 6
9 
10  !> Data structure to hold a global cell identifier,
11  !! using a pointer to the model and its local cell
12  !< index
13  type, public :: globalcelltype
14  integer(I4B) :: index !< the index on the model grid
15  class(virtualmodeltype), pointer :: v_model => null() !< virtual model
16  end type
17 
18  ! a global cell with neighbors
19  type, public :: cellwithnbrstype
20  type(globalcelltype) :: cell
21  integer(I4B) :: nrofnbrs = 0
22  type(cellwithnbrstype), dimension(:), pointer, &
23  contiguous :: neighbors => null()
24  contains
25  procedure :: addnbrcell
26  end type
27 
28 contains
29 
30  subroutine addnbrcell(this, index, v_model)
31  class(cellwithnbrstype) :: this
32  integer(I4B) :: index
33  class(virtualmodeltype), pointer :: v_model
34  ! local
35  integer(I4B) :: nbrCnt, currentSize, i
36  type(cellwithnbrstype), dimension(:), pointer, contiguous :: newNeighbors
37  type(cellwithnbrstype), dimension(:), pointer, contiguous :: oldNeighbors
38 
39  if (.not. associated(this%neighbors)) then
40  allocate (this%neighbors(defaultcapacity))
41  this%nrOfNbrs = 0
42  end if
43 
44  nbrcnt = this%nrOfNbrs
45  currentsize = size(this%neighbors)
46  if (nbrcnt + 1 > currentsize) then
47 
48  ! inflate
49  oldneighbors => this%neighbors
50  allocate (newneighbors(currentsize + defaultcapacity))
51  do i = 1, currentsize
52  newneighbors(i) = oldneighbors(i)
53  end do
54  this%neighbors => newneighbors
55 
56  ! clean up
57  deallocate (oldneighbors)
58  nullify (oldneighbors)
59  end if
60 
61  this%neighbors(nbrcnt + 1)%cell%index = index
62  this%neighbors(nbrcnt + 1)%cell%v_model => v_model
63 
64  this%nrOfNbrs = nbrcnt + 1
65 
66  end subroutine addnbrcell
67 
68 end module
subroutine addnbrcell(this, index, v_model)
integer(i4b), parameter defaultcapacity
index
Definition: CellWithNbrs.f90:8
This module defines variable data types.
Definition: kind.f90:8
Data structure to hold a global cell identifier, using a pointer to the model and its local cell.