MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
sparsemodule Module Reference

Data Types

type  rowtype
 
type  sparsematrix
 

Functions/Subroutines

subroutine destroy (this)
 
subroutine initialize (this, nrow, ncol, rowmaxnnz)
 
subroutine initializefixed (this, nrow, ncol, maxnnz)
 
subroutine filliaja (this, ia, ja, ierror, sort)
 
subroutine addconnection (this, i, j, inodup, iaddop)
 
subroutine insert (j, thisrow, inodup, iadded)
 
subroutine sort (this, with_csr)
 
subroutine sortintarray (nval, iarray)
 
subroutine csr_diagsum (ia, flowja)
 

Function/Subroutine Documentation

◆ addconnection()

subroutine sparsemodule::addconnection ( class(sparsematrix), intent(inout)  this,
integer(i4b), intent(in)  i,
integer(i4b), intent(in)  j,
integer(i4b), intent(in)  inodup,
integer(i4b), intent(inout), optional  iaddop 
)

Definition at line 149 of file Sparse.f90.

150  !add a connection to the sparsematrix. if inodup
151  !(for no duplicates) is 1, then j is added only
152  !if it is unique.
153  ! -- dummy
154  class(sparsematrix), intent(inout) :: this
155  integer(I4B), intent(in) :: i, j, inodup
156  integer(I4B), optional, intent(inout) :: iaddop
157  ! -- local
158  integer(I4B) :: irow_local
159  integer(I4B) :: iadded
160  ! -- code
161  !
162  ! -- when distributed system, reduce row numbers to local range
163  irow_local = i - this%offset
164  !
165  call insert(j, this%row(irow_local), inodup, iadded)
166  this%nnz = this%nnz + iadded
167  if (j < this%offset + 1 .or. j > this%offset + this%nrow) then
168  ! count the off-diagonal entries separately
169  this%nnz_od = this%nnz_od + iadded
170  end if
171  if (present(iaddop)) iaddop = iadded
172  !
173  ! -- return
174  return
Here is the call graph for this function:

◆ csr_diagsum()

subroutine sparsemodule::csr_diagsum ( integer(i4b), dimension(:), contiguous  ia,
real(dp), dimension(:), contiguous  flowja 
)

Definition at line 280 of file Sparse.f90.

281  !Add up the off diagonal terms and put the sum in the
282  !diagonal position
283  ! -- dummy
284  integer(I4B), dimension(:), contiguous :: ia
285  real(DP), dimension(:), contiguous :: flowja
286  ! -- local
287  integer(I4B) :: nodes
288  integer(I4B) :: n
289  integer(I4B) :: iposdiag
290  integer(I4B) :: ipos
291  ! -- code
292  nodes = size(ia) - 1
293  do n = 1, nodes
294  iposdiag = ia(n)
295  do ipos = ia(n) + 1, ia(n + 1) - 1
296  flowja(iposdiag) = flowja(iposdiag) + flowja(ipos)
297  end do
298  end do
299  return
Here is the caller graph for this function:

◆ destroy()

subroutine sparsemodule::destroy ( class(sparsematrix), intent(inout)  this)

Definition at line 34 of file Sparse.f90.

35  class(sparsematrix), intent(inout) :: this
36  deallocate (this%row)

◆ filliaja()

subroutine sparsemodule::filliaja ( class(sparsematrix), intent(inout)  this,
integer(i4b), dimension(:), intent(inout)  ia,
integer(i4b), dimension(:), intent(inout)  ja,
integer(i4b), intent(inout)  ierror,
logical, intent(in), optional  sort 
)

Definition at line 90 of file Sparse.f90.

91  !allocate and fill the ia and ja arrays using information
92  !from the sparsematrix.
93  !ierror is returned as:
94  ! 0 if no error
95  ! 1 if ia is not the correct size
96  ! 2 if ja is not the correct size
97  ! 3 if both ia and ja are not correct size
98  ! -- dummy
99  class(sparsematrix), intent(inout) :: this
100  integer(I4B), dimension(:), intent(inout) :: ia, ja
101  integer(I4B), intent(inout) :: ierror
102  logical, intent(in), optional :: sort
103  ! -- local
104  logical :: sortja
105  integer(I4B) :: i, j, ipos
106  ! -- code
107  !
108  ! -- process optional dummy variables
109  if (present(sort)) then
110  sortja = sort
111  else
112  sortja = .false.
113  end if
114  !
115  ! -- initialize error variable
116  ierror = 0
117  !
118  ! -- check for error conditions
119  if (ubound(ia, dim=1) /= this%nrow + 1) then
120  ierror = 1
121  end if
122  if (ubound(ja, dim=1) /= this%nnz) then
123  ierror = ierror + 2
124  end if
125  if (ierror /= 0) then
126  return
127  end if
128  !
129  ! -- sort this
130  if (sortja) then
131  call this%sort()
132  end if
133  !
134  ! -- fill ia and ja
135  ipos = 1
136  ia(1) = ipos
137  do i = 1, this%nrow
138  do j = 1, this%row(i)%nnz
139  ja(ipos) = this%row(i)%icolarray(j)
140  ipos = ipos + 1
141  end do
142  ia(i + 1) = ipos
143  end do
144  !
145  ! -- return
146  return

◆ initialize()

subroutine sparsemodule::initialize ( class(sparsematrix), intent(inout)  this,
integer(i4b), intent(in)  nrow,
integer(i4b), intent(in)  ncol,
integer(i4b), dimension(nrow), intent(in)  rowmaxnnz 
)

Definition at line 39 of file Sparse.f90.

40  !initial the sparse matrix. This subroutine
41  !acts a method for a sparse matrix by initializing
42  !the row data. It presently requires one maximum
43  !value for all rows, however, this can be changed
44  !to a vector of maximum values with one value for
45  !each row.
46  ! -- dummy
47  class(sparsematrix), intent(inout) :: this
48  integer(I4B), intent(in) :: nrow, ncol
49  integer(I4B), intent(in), dimension(nrow) :: rowmaxnnz
50  ! -- local
51  integer(I4B) :: i
52  ! -- code
53  this%offset = 0
54  this%nrow = nrow
55  this%ncol = ncol
56  this%nnz = 0
57  this%nnz_od = 0
58  allocate (this%row(nrow))
59  do i = 1, nrow
60  allocate (this%row(i)%icolarray(rowmaxnnz(i)))
61  this%row(i)%icolarray = 0
62  this%row(i)%nnz = 0
63  end do
64  !
65  ! -- return
66  return
Here is the caller graph for this function:

◆ initializefixed()

subroutine sparsemodule::initializefixed ( class(sparsematrix), intent(inout)  this,
integer(i4b), intent(in)  nrow,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  maxnnz 
)

Definition at line 70 of file Sparse.f90.

71  implicit none
72  class(sparsematrix), intent(inout) :: this
73  integer(I4B), intent(in) :: nrow, ncol
74  integer(I4B), intent(in) :: maxnnz
75  ! local
76  integer(I4B), dimension(:), allocatable :: rowmaxnnz
77  integer(I4B) :: i
78 
79  allocate (rowmaxnnz(nrow))
80 
81  do i = 1, nrow
82  rowmaxnnz(i) = maxnnz
83  end do
84 
85  call this%initialize(nrow, ncol, rowmaxnnz)
86  deallocate (rowmaxnnz)
87 
Here is the caller graph for this function:

◆ insert()

subroutine sparsemodule::insert ( integer(i4b), intent(in)  j,
type(rowtype), intent(inout)  thisrow,
integer(i4b), intent(in)  inodup,
integer(i4b), intent(inout)  iadded 
)

Definition at line 177 of file Sparse.f90.

178  !insert j into thisrow (for row i)
179  !inodup=1 means do not include duplicate connections
180  !iadded is 1 if a new entry was added (meaning that nnz for the row was increased)
181  !iadded is 0 if duplicate and not added. Used to track total number of connections
182  ! -- dummy
183  integer(I4B), intent(in) :: j, inodup
184  type(rowtype), intent(inout) :: thisrow
185  integer(I4B), allocatable, dimension(:) :: iwk
186  integer(I4B), intent(inout) :: iadded
187  ! -- local
188  integer(I4B) :: jj, maxnnz
189  ! -- code
190  iadded = 0
191  maxnnz = ubound(thisrow%icolarray, dim=1)
192  if (thisrow%icolarray(1) == 0) then
193  thisrow%icolarray(1) = j
194  thisrow%nnz = thisrow%nnz + 1
195  iadded = 1
196  return
197  end if
198  if (thisrow%nnz == maxnnz) then
199  ! -- increase size of the row
200  allocate (iwk(thisrow%nnz))
201  iwk = thisrow%icolarray
202  deallocate (thisrow%icolarray)
203  ! -- Specify how to increase the size of the icolarray. Adding 1
204  ! will be most memory conservative, but may be a little slower
205  ! due to frequent allocate/deallocate. Another option would be
206  ! to double the size: maxnnz=maxnnz*2
207  maxnnz = maxnnz + 1
208  allocate (thisrow%icolarray(maxnnz))
209  thisrow%icolarray(1:thisrow%nnz) = iwk(1:thisrow%nnz)
210  thisrow%icolarray(thisrow%nnz + 1:maxnnz) = 0
211  end if
212  if (inodup == 1) then
213  do jj = 1, thisrow%nnz
214  if (thisrow%icolarray(jj) == j) then
215  return
216  end if
217  end do
218  end if
219  !
220  ! -- add the connection to end
221  thisrow%nnz = thisrow%nnz + 1
222  thisrow%icolarray(thisrow%nnz) = j
223  iadded = 1
224  !
225  ! -- return
226  return
Here is the caller graph for this function:

◆ sort()

subroutine sparsemodule::sort ( class(sparsematrix), intent(inout)  this,
logical(lgp), optional  with_csr 
)

Definition at line 229 of file Sparse.f90.

230  !sort the icolarray for each row, but do not include
231  !the diagonal position in the sort so that it stays in front
232  ! -- dummy
233  class(sparsematrix), intent(inout) :: this
234  logical(LGP), optional :: with_csr
235  ! -- local
236  integer(I4B) :: i, nval, start_idx
237  ! -- code
238  start_idx = 2
239  if (present(with_csr)) then
240  if (with_csr) then
241  ! CSR: don't put diagonal up front
242  start_idx = 1
243  end if
244  end if
245 
246  do i = 1, this%nrow
247  nval = this%row(i)%nnz
248  call sortintarray(nval - start_idx + 1, &
249  this%row(i)%icolarray(start_idx:nval))
250  end do
251  !
252  ! -- return
253  return
Here is the call graph for this function:

◆ sortintarray()

subroutine sparsemodule::sortintarray ( integer(i4b), intent(in)  nval,
integer(i4b), dimension(nval), intent(inout)  iarray 
)

Definition at line 256 of file Sparse.f90.

257  !simple subroutine for sorting an array
258  !in place. It is not the fastest sort function
259  !but should suffice for relatively short nodelists.
260  ! -- dummy
261  integer(I4B), intent(in) :: nval
262  integer(I4B), intent(inout), dimension(nval) :: iarray
263  ! -- local
264  integer(I4B) :: i, j, itemp
265  ! -- code
266  do i = 1, nval - 1
267  do j = i + 1, nval
268  if (iarray(i) > iarray(j)) then
269  itemp = iarray(j)
270  iarray(j) = iarray(i)
271  iarray(i) = itemp
272  end if
273  end do
274  end do
275  !
276  ! -- return
277  return
Here is the caller graph for this function: