MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
sparsematrixmodule Module Reference

Data Types

type  sparsematrixtype
 

Functions/Subroutines

subroutine spm_init (this, sparse, mem_path)
 Initialize sparse matrix from passed. More...
 
subroutine spm_destroy (this)
 
class(vectorbasetype) function, pointer spm_create_vec_mm (this, n, name, mem_path)
 
class(vectorbasetype) function, pointer spm_create_vec (this, n)
 
real(dp) function spm_get_value_pos (this, ipos)
 
real(dp) function spm_get_diag_value (this, irow)
 
subroutine spm_set_diag_value (this, irow, diag_value)
 
subroutine spm_set_value_pos (this, ipos, value)
 
subroutine spm_add_value_pos (this, ipos, value)
 
subroutine spm_add_diag_value (this, irow, value)
 
integer(i4b) function spm_get_first_col_pos (this, irow)
 
integer(i4b) function spm_get_last_col_pos (this, irow)
 
integer(i4b) function spm_get_column (this, ipos)
 
integer(i4b) function spm_get_position (this, irow, icol)
 Return position index for (irow,icol) element in the matrix. This index can be used in other routines for direct access. More...
 
integer(i4b) function spm_get_position_diag (this, irow)
 
subroutine allocate_scalars (this)
 
subroutine allocate_arrays (this)
 
subroutine spm_zero_entries (this)
 Set all entries in the matrix to zero. More...
 
subroutine spm_zero_row_offdiag (this, irow)
 Set all off-diagonal entries in the matrix to zero. More...
 
subroutine spm_get_aij (this, ia, ja, amat)
 
integer(i4b) function spm_get_row_offset (this)
 
subroutine spm_multiply (this, vec_x, vec_y)
 Calculates the matrix vector product y = A*x. More...
 

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine sparsematrixmodule::allocate_arrays ( class(sparsematrixtype this)
private

Definition at line 240 of file SparseMatrix.f90.

241  class(SparseMatrixType) :: this
242 
243  call mem_allocate(this%ia, this%nrow + 1, 'IA', this%memory_path)
244  call mem_allocate(this%ja, this%nja, 'JA', this%memory_path)
245  call mem_allocate(this%amat, this%nja, 'AMAT', this%memory_path)
246 

◆ allocate_scalars()

subroutine sparsematrixmodule::allocate_scalars ( class(sparsematrixtype this)
private

Definition at line 231 of file SparseMatrix.f90.

232  class(SparseMatrixType) :: this
233 
234  call mem_allocate(this%nrow, 'NROW', this%memory_path)
235  call mem_allocate(this%ncol, 'NCOL', this%memory_path)
236  call mem_allocate(this%nja, 'NJA', this%memory_path)
237 

◆ spm_add_diag_value()

subroutine sparsematrixmodule::spm_add_diag_value ( class(sparsematrixtype this,
integer(i4b)  irow,
real(dp)  value 
)
private

Definition at line 162 of file SparseMatrix.f90.

163  class(SparseMatrixType) :: this
164  integer(I4B) :: irow
165  real(DP) :: value
166 
167  this%amat(this%ia(irow)) = this%amat(this%ia(irow)) + value
168 

◆ spm_add_value_pos()

subroutine sparsematrixmodule::spm_add_value_pos ( class(sparsematrixtype this,
integer(i4b)  ipos,
real(dp)  value 
)
private

Definition at line 153 of file SparseMatrix.f90.

154  class(SparseMatrixType) :: this
155  integer(I4B) :: ipos
156  real(DP) :: value
157 
158  this%amat(ipos) = this%amat(ipos) + value
159 

◆ spm_create_vec()

class(vectorbasetype) function, pointer sparsematrixmodule::spm_create_vec ( class(sparsematrixtype this,
integer(i4b)  n 
)
private
Parameters
nthe nr. of elements in the vector

Definition at line 104 of file SparseMatrix.f90.

105  class(SparseMatrixType) :: this ! this sparse matrix
106  integer(I4B) :: n !< the nr. of elements in the vector
107  class(VectorBaseType), pointer :: vec ! the vector to create
108  ! local
109  class(SeqVectorType), pointer :: seq_vec
110 
111  allocate (seq_vec)
112  call seq_vec%create(n)
113  vec => seq_vec
114 

◆ spm_create_vec_mm()

class(vectorbasetype) function, pointer sparsematrixmodule::spm_create_vec_mm ( class(sparsematrixtype this,
integer(i4b)  n,
character(len=*)  name,
character(len=*)  mem_path 
)
private
Parameters
nthe nr. of elements in the vector
namethe variable name (for access through memory manager)
mem_pathmemory path for storing the underlying memory items

Definition at line 89 of file SparseMatrix.f90.

90  class(SparseMatrixType) :: this ! this sparse matrix
91  integer(I4B) :: n !< the nr. of elements in the vector
92  character(len=*) :: name !< the variable name (for access through memory manager)
93  character(len=*) :: mem_path !< memory path for storing the underlying memory items
94  class(VectorBaseType), pointer :: vec ! the vector to create
95  ! local
96  class(SeqVectorType), pointer :: seq_vec
97 
98  allocate (seq_vec)
99  call seq_vec%create_mm(n, name, mem_path)
100  vec => seq_vec
101 

◆ spm_destroy()

subroutine sparsematrixmodule::spm_destroy ( class(sparsematrixtype this)
private

Definition at line 76 of file SparseMatrix.f90.

77  class(SparseMatrixType) :: this
78 
79  call mem_deallocate(this%nrow)
80  call mem_deallocate(this%ncol)
81  call mem_deallocate(this%nja)
82 
83  call mem_deallocate(this%ia)
84  call mem_deallocate(this%ja)
85  call mem_deallocate(this%amat)
86 

◆ spm_get_aij()

subroutine sparsematrixmodule::spm_get_aij ( class(sparsematrixtype this,
integer(i4b), dimension(:), pointer, contiguous  ia,
integer(i4b), dimension(:), pointer, contiguous  ja,
real(dp), dimension(:), pointer, contiguous  amat 
)
private

Definition at line 276 of file SparseMatrix.f90.

277  class(SparseMatrixType) :: this
278  integer(I4B), dimension(:), pointer, contiguous :: ia
279  integer(I4B), dimension(:), pointer, contiguous :: ja
280  real(DP), dimension(:), pointer, contiguous :: amat
281 
282  ia => this%ia
283  ja => this%ja
284  amat => this%amat
285 

◆ spm_get_column()

integer(i4b) function sparsematrixmodule::spm_get_column ( class(sparsematrixtype this,
integer(i4b)  ipos 
)
private

Definition at line 189 of file SparseMatrix.f90.

190  class(SparseMatrixType) :: this
191  integer(I4B) :: ipos
192  integer(I4B) :: icol
193 
194  icol = this%ja(ipos)
195 

◆ spm_get_diag_value()

real(dp) function sparsematrixmodule::spm_get_diag_value ( class(sparsematrixtype this,
integer(i4b)  irow 
)
private

Definition at line 126 of file SparseMatrix.f90.

127  class(SparseMatrixType) :: this
128  integer(I4B) :: irow
129  real(DP) :: diag_value
130 
131  diag_value = this%amat(this%ia(irow))
132 

◆ spm_get_first_col_pos()

integer(i4b) function sparsematrixmodule::spm_get_first_col_pos ( class(sparsematrixtype this,
integer(i4b)  irow 
)
private

Definition at line 171 of file SparseMatrix.f90.

172  class(SparseMatrixType) :: this
173  integer(I4B) :: irow
174  integer(I4B) :: first_col_pos
175 
176  first_col_pos = this%ia(irow)
177 

◆ spm_get_last_col_pos()

integer(i4b) function sparsematrixmodule::spm_get_last_col_pos ( class(sparsematrixtype this,
integer(i4b)  irow 
)
private

Definition at line 180 of file SparseMatrix.f90.

181  class(SparseMatrixType) :: this
182  integer(I4B) :: irow
183  integer(I4B) :: last_col_pos
184 
185  last_col_pos = this%ia(irow + 1) - 1
186 

◆ spm_get_position()

integer(i4b) function sparsematrixmodule::spm_get_position ( class(sparsematrixtype this,
integer(i4b)  irow,
integer(i4b)  icol 
)
private

Definition at line 202 of file SparseMatrix.f90.

203  class(SparseMatrixType) :: this
204  integer(I4B) :: irow
205  integer(I4B) :: icol
206  integer(I4B) :: ipos
207  ! local
208  integer(I4B) :: i, icol_s, icol_e
209 
210  ipos = -1
211  icol_s = this%get_first_col_pos(irow)
212  icol_e = this%get_last_col_pos(irow)
213  do i = icol_s, icol_e
214  if (this%ja(i) == icol) then
215  ipos = i
216  return
217  end if
218  end do
219 

◆ spm_get_position_diag()

integer(i4b) function sparsematrixmodule::spm_get_position_diag ( class(sparsematrixtype this,
integer(i4b)  irow 
)
private

Definition at line 222 of file SparseMatrix.f90.

223  class(SparseMatrixType) :: this
224  integer(I4B) :: irow
225  integer(I4B) :: ipos_diag
226 
227  ipos_diag = this%ia(irow)
228 

◆ spm_get_row_offset()

integer(i4b) function sparsematrixmodule::spm_get_row_offset ( class(sparsematrixtype this)
private

Definition at line 288 of file SparseMatrix.f90.

289  class(SparseMatrixType) :: this
290  integer(I4B) :: offset
291 
292  offset = 0
293 

◆ spm_get_value_pos()

real(dp) function sparsematrixmodule::spm_get_value_pos ( class(sparsematrixtype this,
integer(i4b)  ipos 
)
private

Definition at line 117 of file SparseMatrix.f90.

118  class(SparseMatrixType) :: this
119  integer(I4B) :: ipos
120  real(DP) :: value
121 
122  value = this%amat(ipos)
123 

◆ spm_init()

subroutine sparsematrixmodule::spm_init ( class(sparsematrixtype this,
type(sparsematrix sparse,
character(len=*)  mem_path 
)

Definition at line 54 of file SparseMatrix.f90.

55  class(SparseMatrixType) :: this
56  type(sparsematrix) :: sparse
57  character(len=*) :: mem_path
58  ! local
59  integer(I4B) :: ierror
60 
61  this%memory_path = mem_path
62 
63  call this%allocate_scalars()
64 
65  this%nrow = sparse%nrow
66  this%ncol = sparse%ncol
67  this%nja = sparse%nnz
68 
69  call this%allocate_arrays()
70 
71  call sparse%filliaja(this%ia, this%ja, ierror, sort=.false.)
72  call this%zero_entries()
73 

◆ spm_multiply()

subroutine sparsematrixmodule::spm_multiply ( class(sparsematrixtype this,
class(vectorbasetype), pointer  vec_x,
class(vectorbasetype), pointer  vec_y 
)
private

Definition at line 298 of file SparseMatrix.f90.

299  class(SparseMatrixType) :: this
300  class(VectorBaseType), pointer :: vec_x
301  class(VectorBaseType), pointer :: vec_y
302  ! local
303  integer(I4B) :: irow, icol, ipos
304  real(DP), dimension(:), pointer, contiguous :: x, y
305 
306  x => vec_x%get_array()
307  y => vec_y%get_array()
308  do irow = 1, this%nrow
309  y(irow) = dzero
310  do ipos = this%ia(irow), this%ia(irow + 1) - 1
311  icol = this%ja(ipos)
312  y(irow) = y(irow) + this%amat(ipos) * x(icol)
313  end do
314  end do
315 

◆ spm_set_diag_value()

subroutine sparsematrixmodule::spm_set_diag_value ( class(sparsematrixtype this,
integer(i4b)  irow,
real(dp)  diag_value 
)
private

Definition at line 135 of file SparseMatrix.f90.

136  class(SparseMatrixType) :: this
137  integer(I4B) :: irow
138  real(DP) :: diag_value
139 
140  this%amat(this%ia(irow)) = diag_value
141 

◆ spm_set_value_pos()

subroutine sparsematrixmodule::spm_set_value_pos ( class(sparsematrixtype this,
integer(i4b)  ipos,
real(dp)  value 
)
private

Definition at line 144 of file SparseMatrix.f90.

145  class(SparseMatrixType) :: this
146  integer(I4B) :: ipos
147  real(DP) :: value
148 
149  this%amat(ipos) = value
150 

◆ spm_zero_entries()

subroutine sparsematrixmodule::spm_zero_entries ( class(sparsematrixtype this)
private

Definition at line 251 of file SparseMatrix.f90.

252  class(SparseMatrixType) :: this
253  ! local
254  integer(I4B) :: i
255 
256  do i = 1, this%nja
257  this%amat(i) = dzero
258  end do
259 

◆ spm_zero_row_offdiag()

subroutine sparsematrixmodule::spm_zero_row_offdiag ( class(sparsematrixtype this,
integer(i4b)  irow 
)
private

Definition at line 264 of file SparseMatrix.f90.

265  class(SparseMatrixType) :: this
266  integer(I4B) :: irow
267  ! local
268  integer(I4B) :: ipos
269 
270  do ipos = this%ia(irow) + 1, this%ia(irow + 1) - 1
271  this%amat(ipos) = dzero
272  end do
273