MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
BaseMatrix_module.f90
Go to the documentation of this file.
1! Copyright 2019
2!
3! For a comprehensive list of the developers that contributed to these codes
4! see the UK-AMOR website.
5!
6! This file is part of UKRmol-in (UKRmol+ suite).
7!
8! UKRmol-in is free software: you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation, either version 3 of the License, or
11! (at your option) any later version.
12!
13! UKRmol-in is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17!
18! You should have received a copy of the GNU General Public License
19! along with UKRmol-in (in source/COPYING). Alternatively, you can also visit
20! <https://www.gnu.org/licenses/>.
21
22!> \brief Base matrix module
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
27!>
28module basematrix_module
29
30 use const_gbl, only: stdout
32 use options_module, only: options
33 use precisn, only: wp
34
35 implicit none
36
37 public basematrix
38
39 private
40
41 !> \brief Base matrix type
42 !> \authors A Al-Refaie
43 !> \date 2017
44 !>
45 !> This handles the matrix elements and also expands the vector size if we have reached max capacity.
46 !>
47 type, abstract :: basematrix
48 real(wp) :: threshold
49 integer :: matrix_type
50 integer :: matrix_indexing = matrix_index_fortran
51 integer :: matrix_mode
52 integer :: excluded_rowcol = 2000000000
53 !>Number of elements in the array of both the integral and coefficients
54 logical :: constructed = .false.
55 integer :: matrix_dimension
56 integer :: n = 0
57 integer :: block_size
58 logical :: remove_row_column = .false.
59 real(wp), allocatable :: diagonal(:)
60 logical :: initialized = .false.
61 contains
62 procedure, public :: initialize_matrix_structure
63 procedure, public :: construct
64 procedure, public :: insert_matrix_element
65 procedure, public :: get_matrix_element
66 procedure, public :: exclude_row_column
67 procedure, public :: is_empty
68 procedure, public :: get_size
69 procedure, public :: get_matrix_size
70
71 procedure, public :: clear
72 procedure, public :: destroy
73
74 procedure, public :: update_continuum
75 procedure, public :: set_options
76 procedure, public :: update_pure_l2
77 procedure, public :: finalize_matrix
78 procedure, public :: store_diagonal
79 procedure, public :: initialize_struct_self
80
81 procedure(generic_construct), deferred :: construct_self
82 procedure(generic_insert), deferred :: insert_matelem_self
83 procedure(generic_get), deferred :: get_matelem_self
84 procedure(generic_clear), deferred :: clear_self
85 procedure(generic_destroy), deferred :: destroy_self
86
87 procedure, public :: expand_capacity
88 procedure, public :: print => print_mat
89
90 procedure, private :: check_bounds
91 procedure, private :: check_constructed
92 end type
93
94 abstract interface
95 subroutine generic_construct (this)
96 import :: basematrix
97 class(basematrix) :: this
98 end subroutine generic_construct
99 end interface
100
101 abstract interface
102 subroutine generic_insert(this, i, j, coefficient, class, thresh)
103 use precisn, only: wp
104 import :: basematrix
105 class(basematrix) :: this
106 integer, intent(in) :: i, j
107 real(wp), intent(in) :: coefficient
108 integer, intent(in) :: class
109 real(wp), intent(in) :: thresh
110 end subroutine generic_insert
111 end interface
112
113 abstract interface
114 subroutine generic_get (this, idx, i, j, coeff)
115 use precisn, only: wp
116 import :: basematrix
117 class(basematrix) :: this
118 integer, intent(in) :: idx
119 integer, intent(out) :: i, j
120 real(wp), intent(out) :: coeff
121 end subroutine generic_get
122 end interface
123
124 abstract interface
125 subroutine generic_clear (this)
126 import :: basematrix
127 class(basematrix) :: this
128 end subroutine generic_clear
129 end interface
130
131 abstract interface
132 subroutine generic_destroy (this)
133 import :: basematrix
134 class(basematrix) :: this
135 end subroutine generic_destroy
136 end interface
137
138contains
139
140 !> \brief Check that matrix is constructed
141 !> \authors A Al-Refaie
142 !> \date 2017
143 !>
144 !> Check that matrix is constructed; hard stop if not.
145 !>
146 subroutine check_constructed (this)
147 class(basematrix) :: this
148
149 if (.not. this % constructed) then
150 write (stdout, "('Vector::constructed - Vector is not constructed')")
151 stop "Vector::constructed - Vector is not constructed"
152 end if
153
154 end subroutine check_constructed
155
156
157 !> \brief Construct the matrix
158 !> \authors A Al-Refaie
159 !> \date 2017
160 !>
161 subroutine construct (this, mat_mode, threshold)
162 class(basematrix) :: this
163 integer, intent(in) :: mat_mode
164 real(wp), optional, intent(in) :: threshold
165
166 integer :: err
167
168 if (present(threshold)) then
169 this % threshold = threshold
170 else
171 this % threshold = default_matelem_threshold
172 end if
173
174 this % matrix_mode = mat_mode
175
176 call this % construct_self
177 this % remove_row_column = .false.
178 this % constructed = .true.
179
180 end subroutine construct
181
182
183 !> \brief ?
184 !> \authors A Al-Refaie
185 !> \date 2017
186 !>
187 subroutine exclude_row_column (this, row_column)
188 class(basematrix) :: this
189 integer, intent(in) :: row_column
190
191 this % excluded_rowcol = row_column
192 if (this % excluded_rowcol > 0) then
193 this % remove_row_column = .true.
194 end if
195
196 end subroutine exclude_row_column
197
198
199 !> \brief ?
200 !> \authors A Al-Refaie
201 !> \date 2017
202 !>
203 subroutine initialize_matrix_structure (this, matrix_size, matrix_type, block_size)
204 class(basematrix) :: this
205 integer, intent(in) :: matrix_size, matrix_type, block_size
206 integer :: ido, ierr
207
208 this % initialized = .true.
209 this % matrix_dimension = matrix_size
210 this % matrix_type = matrix_type
211 this % block_size = block_size
212
213 if (this % remove_row_column .and. this % excluded_rowcol <= this % matrix_dimension) then
214 this % matrix_dimension = this % matrix_dimension - 1
215 end if
216
217 if (allocated(this % diagonal)) deallocate(this % diagonal)
218 write (stdout, "('Allocating dimension ',i12)") this % matrix_dimension
219
220 allocate(this % diagonal(this % matrix_dimension), stat = ierr)
221 this % diagonal = 0
222 if (ierr /= 0) then
223 stop "ERROR"
224 end if
225
226 call this % initialize_struct_self(this % matrix_dimension, this % matrix_type, this % block_size)
227
228 end subroutine initialize_matrix_structure
229
230
231 !> \brief ?
232 !> \authors A Al-Refaie
233 !> \date 2017
234 !>
235 !> Usually does nothing but some matrix formats may need more parameters in order to function.
236 !>
237 subroutine set_options (this, option_val)
238 class(basematrix) :: this
239 class(options), intent(in) :: option_val
240 end subroutine set_options
241
242
243 !> \brief Initialize the type
244 !> \authors A Al-Refaie
245 !> \date 2017
246 !>
247 !> Do nothing
248 !>
249 subroutine initialize_struct_self (this, matrix_size, matrix_type, block_size)
250 class(basematrix) :: this
251 integer, intent(in) :: matrix_size, matrix_type, block_size
252 end subroutine initialize_struct_self
253
254
255 !> \brief Set matrix element
256 !> \authors A Al-Refaie
257 !> \date 2017
258 !>
259 subroutine insert_matrix_element (this, i, j, coefficient, class_, thresh_)
260 class(basematrix) :: this
261 integer, intent(in) :: i, j
262 real(wp), intent(in) :: coefficient
263 integer, intent(in), optional :: class_
264 real(wp), intent(in), optional :: thresh_
265 real(wp) :: thresh
266 integer :: class, i_, j_
267
268 !this%matrix_size = max(this%matrix_size,i,j)
269
270 if (.not. this % initialized) then
271 stop "Matrix structure not initialized"
272 end if
273
274 if (present(thresh_)) then
275 thresh = thresh_
276 else
277 thresh = this % threshold
278 endif
279
280 if (present(class_)) then
281 class = class_
282 else
283 class = -1
284 end if
285
286 !Make sure we are in the lower triangular
287 if (i >= j) then
288 i_ = i
289 j_ = j
290 else
291 i_ = j
292 j_ = i
293 end if
294
295 !Exclude if needed
296 if (this % remove_row_column) then
297 if (i_ == this % excluded_rowcol) return
298 if (j_ == this % excluded_rowcol) return
299 if (i_ > this % excluded_rowcol) i_ = i_ - 1
300 if (j_ > this % excluded_rowcol) j_ = j_ - 1
301 end if
302 if (i == j) call this % store_diagonal(i, coefficient)
303 if (abs(coefficient) < thresh) return
304
305 call this % insert_matelem_self(i_, j_, coefficient, class, thresh)
306
307 end subroutine insert_matrix_element
308
309
310 !> \brief Set diagonal element
311 !> \authors A Al-Refaie
312 !> \date 2017
313 !>
314 subroutine store_diagonal (this, i, coeff)
315 class(basematrix) :: this
316 integer, intent(in) :: i
317 real(wp), intent(in) :: coeff
318
319 !if(myrank /= master) return
320 this % diagonal(i) = coeff
321
322 end subroutine store_diagonal
323
324
325 !> \brief Retrieve matrix element
326 !> \authors A Al-Refaie
327 !> \date 2017
328 !>
329 subroutine get_matrix_element (this, idx, i, j, coeff)
330 class(basematrix) :: this
331 integer, intent(in) :: idx
332 integer, intent(out) :: i, j
333 real(wp), intent(out) :: coeff
334
335 if (this % check_bounds(idx)) then
336 call this % get_matelem_self(idx, i, j, coeff)
337 end if
338
339 end subroutine get_matrix_element
340
341
342 !> \brief ?
343 !> \authors A Al-Refaie
344 !> \date 2017
345 !>
346 subroutine update_continuum (this, force_update)
347 class(basematrix) :: this
348 logical, intent(in) :: force_update
349 end subroutine update_continuum
350
351
352 !> \brief ?
353 !> \authors A Al-Refaie
354 !> \date 2017
355 !>
356 subroutine update_pure_l2 (this, force_update, count_)
357 class(basematrix) :: this
358 logical, intent(in) :: force_update
359 integer, optional :: count_
360 end subroutine update_pure_l2
361
362
363 !> \brief ?
364 !> \authors A Al-Refaie
365 !> \date 2017
366 !>
367 subroutine finalize_matrix (this)
368 class(basematrix) :: this
369 end subroutine finalize_matrix
370
371
372 !> \brief ?
373 !> \authors A Al-Refaie
374 !> \date 2017
375 !>
376 subroutine expand_capacity (this, capacity)
377 class(basematrix) :: this
378 integer, intent(in) :: capacity
379 end subroutine expand_capacity
380
381
382 !> \brief ?
383 !> \authors A Al-Refaie
384 !> \date 2017
385 !>
386 logical function check_bounds (this, i)
387 class(basematrix) :: this
388 integer, intent(in) :: i
389
390 if (i <= 0 .or. i > this % n) then
391 write (stdout, "('Vector::check_bounds - Out of Bounds access')")
392 stop "Vector::check_bounds - Out of Bounds access"
393 check_bounds = .false.
394 else
395 check_bounds = .true.
396 end if
397
398 end function check_bounds
399
400
401 !> \brief Determine if matrix is empty
402 !> \authors A Al-Refaie
403 !> \date 2017
404 !>
405 !> Check if there are no elements in the mastrix.
406 !>
407 logical function is_empty(this)
408 class(basematrix) :: this
409
410 is_empty = (this % n == 0)
411
412 end function is_empty
413
414
415 !> \brief Clear matrix
416 !> \authors A Al-Refaie
417 !> \date 2017
418 !>
419 subroutine clear (this)
420 class(basematrix) :: this
421 integer :: ido
422
423 !do ido=1,this%max_capacity
424 ! this%matrix_elements(ido)%i = -1
425 ! this%matrix_elements(ido)%j = -1
426 ! this%matrix_elements(ido)%coefficient = 0.0_wp
427 !enddo
428 !this%matrix_size = 0
429
430 this % n = 0
431 call this % clear_self
432 this % initialized = .false.
433
434 end subroutine clear
435
436
437 !> \brief Get matrix size (rank)
438 !> \authors A Al-Refaie
439 !> \date 2017
440 !>
441 integer function get_size (this)
442 class(basematrix) :: this
443
444 get_size = this % n
445
446 end function get_size
447
448
449 !> \brief Get matrix size (number of elements)
450 !> \authors A Al-Refaie
451 !> \date 2017
452 !>
453 integer function get_matrix_size(this)
454 class(basematrix) :: this
455
456 if (.not. this % initialized) then
457 stop "Matrix structure not initialized"
458 end if
459
460 get_matrix_size = this % matrix_dimension
461
462 end function get_matrix_size
463
464
465 !> \brief Print matrix
466 !> \authors A Al-Refaie
467 !> \date 2017
468 !>
469 subroutine print_mat (this)
470 class(basematrix) :: this
471 integer :: i, j, ido
472 real(wp) :: coeff
473
474 do ido = 1, this % n
475 call this % get_matrix_element(ido, i, j, coeff)
476 write (stdout, "(i8,i8,D16.8)") i, j, coeff
477 end do
478
479 end subroutine print_mat
480
481
482 !> \brief Destroy matrix
483 !> \authors A Al-Refaie
484 !> \date 2017
485 !>
486 subroutine destroy (this)
487 class(basematrix) :: this
488 integer :: num_arrays,ido
489
490 if (allocated(this % diagonal)) deallocate(this % diagonal)
491
492 call this % destroy_self
493
494 this % constructed = .false.
495
496 end subroutine destroy
497
498
499 !> \brief Update matrix
500 !> \authors A Al-Refaie
501 !> \date 2017
502 !>
503 subroutine update_matrix (this)
504 class(basematrix) :: this
505 end subroutine update_matrix
506
507end module basematrix_module
MPI SCATCI Constants module.
integer, parameter matrix_index_fortran
Matrix uses FORTRAN indexing.
real(wp), parameter default_matelem_threshold
Matrix element storage threshold.