MPI-SCATCI  2.0
An MPI version of SCATCI
MatrixCache_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 
29 
30  use precisn, only: wp
31  use const_gbl, only: stdout
32  use timing_module, only: master_timer
35 
36  implicit none
37 
38  public matrixcache
39 
40  private
41 
42  type :: matrixarray
43  integer, pointer :: ij(:,:)
44  real(wp), pointer :: coefficient(:)
45  integer :: num_elems
46  logical :: constructed = .false.
47  contains
48  !procedure, public :: is_valid
49  procedure, public :: destroy => destroy_matrix_array
50  procedure, public :: get => get_from_array
51  procedure, public :: insert => insert_into_array
52  procedure, public :: construct => construct_array
53  procedure, public :: clear => clear_array
54  end type matrixarray
55 
60  type :: matrixcache
61  type(matrixarray), allocatable :: matrix_arrays(:)
62  integer :: n
63  logical :: constructed = .false.
64  integer :: matrix_index
65  integer :: num_matrix_chunks = 1
66  integer :: max_capacity = 0
67  integer :: expand_size = default_expand_size
68  contains
69  procedure, public :: construct
70  procedure, public :: insert_into_cache
71  procedure, public :: get_from_cache
72  procedure, public :: is_empty
73  procedure, public :: clear
74  procedure, public :: get_size
75  procedure, public :: expand_capacity
76  !procedure, public :: prune_threshold
77  procedure, public :: sort => qsort_matelem
78  procedure, public :: print => print_matelems
79  procedure, public :: get_local_mat
80  procedure, public :: get_chunk_idx
81  procedure, private :: quick_sort_function
82  procedure, private :: partition_function
83  procedure, private :: expand_array
84  procedure, public :: shrink_capacity
85  procedure, public :: clear_and_shrink
86  procedure, private :: check_bounds
87  procedure, public :: destroy => destroy_cache
88  procedure, private :: check_constructed
89  !procedure, public :: count_occurance
90  !procedure, public :: construct
91  end type matrixcache
92 
93 contains
94 
95  subroutine check_constructed (this)
96  class(matrixcache) :: this
97 
98  if (.not. this % constructed) then
99  write (stdout, "('Vector::constructed - Vector is not constructed')")
100  stop "Vector::constructed - Vector is not constructed"
101  end if
102 
103  end subroutine check_constructed
104 
105 
106  integer function get_chunk_idx (this, idx)
107  class(matrixcache) :: this
108  integer, intent(in) :: idx
109  integer :: mat_idx, local_idx
110 
111  call this % get_local_mat(idx, mat_idx, local_idx)
112 
113  get_chunk_idx = mat_idx
114 
115  end function get_chunk_idx
116 
117 
118  subroutine construct (this, expand_size_)
119  class(matrixcache) :: this
120  integer :: err
121  integer, optional, intent(in) :: expand_size_
122 
123  if (present(expand_size_)) then
124  this % expand_size = expand_size_
125  else
126  this % expand_size = default_expand_size
127  end if
128 
129  this % max_capacity = this % expand_size
130  this % n = 0
131  this % matrix_index = 1
132 
133  !Allocate the vectors
134  allocate(this%matrix_arrays(this%matrix_index),stat=err)
135  call master_memory%track_memory(36,size(this%matrix_arrays),err,'MATRIXCACHE::MATRIXARRAY')
136  !allocate(this%matrix_arrays(this%matrix_index)%ij(2,this%expand_size), &
137  ! & this%matrix_arrays(this%matrix_index)%coefficient(this%expand_size),stat=err)
138  !call master_memory%track_memory(kind(this%matrix_arrays(this%matrix_index)%ij),size(this%matrix_arrays(this%matrix_index)%ij),err,'MATRIXCACHE::MATRIXARRAY::IJ')
139  !call master_memory%track_memory(kind(this%matrix_arrays(this%matrix_index)%coefficient),size(this%matrix_arrays(this%matrix_index)%coefficient),err,'MATRIXCACHE::MATRIXARRAY::COEFF')
140 
141  call this % matrix_arrays(this % matrix_index) % construct(this % expand_size)
142  !this%matrix_arrays(this%matrix_index)%num_elems = 0
143 
144  this % num_matrix_chunks = 1
145 
146  call this % clear
147 
148  if (err /= 0) then
149  write (stdout, "('Matrix Cache::construct- arrays not allocated')")
150  stop "Matrix Cache:: arrays not allocated"
151  end if
152 
153  call this % clear
154  this % constructed = .true.
155 
156  end subroutine construct
157 
158 
159  subroutine insert_into_cache (this, i, j, coefficient)
160  class(matrixcache) :: this
161  integer, intent(in) :: i, j
162  real(wp), intent(in) :: coefficient
163  integer :: mat_id, local_index
164 
165  !this % matrix_size = max(this % matrix_size, i, j)
166  this % n = this % n + 1
167 
168  if (this % n > this % max_capacity) then
169  call this % expand_array()
170  end if
171 
172  call this % get_local_mat(this % n, mat_id, local_index)
173 
174  this % num_matrix_chunks = mat_id
175  this % matrix_arrays(mat_id) % ij(1, local_index) = i
176  this % matrix_arrays(mat_id) % ij(2, local_index) = j
177  this % matrix_arrays(mat_id) % coefficient(local_index) = coefficient
178 
179  this % matrix_arrays(mat_id) % num_elems = this % matrix_arrays(mat_id) % num_elems + 1
180 
181  end subroutine insert_into_cache
182 
183 
184  subroutine get_local_mat (this, idx, mat_id, local_index)
185  class(matrixcache) :: this
186  integer, intent(in) :: idx
187  integer, intent(out) :: mat_id
188  integer, intent(out) :: local_index
189 
190  if (this % check_bounds(idx)) then
191  mat_id = (idx - 1) / this % expand_size + 1
192  !if(mat_id > this%num_matrix_chunks) then
193  !write(stdout,"('Error, referencing a chunk larger than me! [idx,max_Capacity,expand_size,num_chunks,referenced_chunk] = ',5i8)") idx,this%max_capacity,this%expand_size,this%num_matrix_chunks,mat_id
194  !stop "matid error"
195  !endif
196  local_index = idx - (mat_id - 1) * this % expand_size
197  end if
198 
199  end subroutine get_local_mat
200 
201 
202  subroutine get_from_cache (this, idx, i, j, coeff)
203  class(matrixcache) :: this
204  integer, intent(in) :: idx
205  integer, intent(out) :: i, j
206  real(wp), intent(out) :: coeff
207 
208  integer :: mat_id, local_index
209 
210  if (this % check_bounds(idx)) then
211  call this % get_local_mat(idx, mat_id, local_index)
212 
213  i = this % matrix_arrays(mat_id) % ij(1, local_index)
214  j = this % matrix_arrays(mat_id) % ij(2, local_index)
215 
216  coeff = this % matrix_arrays(mat_id) % coefficient(local_index)
217  end if
218 
219  end subroutine get_from_cache
220 
221 
222  subroutine expand_array (this)
223  class(matrixcache) :: this
224  type(matrixarray) :: temp_matrix(this % num_matrix_chunks + 1)
225  integer :: new_num_mats, err
226 
227  call this % check_constructed
228 
229  !call master_timer%start_timer("Expand array")
230  temp_matrix(1:this % num_matrix_chunks) = this % matrix_arrays(1:this % num_matrix_chunks)
231  call master_memory % free_memory(36, size(this % matrix_arrays))
232  deallocate(this % matrix_arrays)
233 
234  this % max_capacity = this % max_capacity + this % expand_size
235 
236  new_num_mats = this % max_capacity / this % expand_size
237  this % matrix_index = new_num_mats
238  this % num_matrix_chunks = this % num_matrix_chunks + 1
239  allocate(this % matrix_arrays(this % num_matrix_chunks), stat = err)
240  call master_memory % track_memory(36, size(this % matrix_arrays), err, 'MATRIXCACHE::MATRIXARRAY')
241 
242  this % matrix_arrays(:) = temp_matrix(:)
243  this % matrix_arrays(this % num_matrix_chunks) % num_elems = 0
244  call this % matrix_arrays(this % num_matrix_chunks) % construct(this % expand_size)
245  !call master_timer%stop_timer("Expand array")
246 
247  end subroutine expand_array
248 
249 
250  subroutine destroy_matrix_array (this)
251  class(matrixarray) :: this
252 
253  if (.not. this % constructed) return
254  if (associated(this % ij)) then
255  call master_memory % free_memory(kind(this % ij), size(this % ij))
256  deallocate(this % ij)
257  this % ij => null()
258  end if
259  if (associated(this % coefficient)) then
260  call master_memory % free_memory(kind(this % coefficient), size(this % coefficient))
261  deallocate(this % coefficient)
262  this % coefficient => null()
263  end if
264  this % constructed = .false.
265  this % num_elems = 0
266 
267  end subroutine destroy_matrix_array
268 
269 
270  subroutine construct_array (this, capacity)
271  class(matrixarray) :: this
272  integer, intent(in) :: capacity
273  integer :: err
274 
275  if (.not. this % constructed) then
276  call master_memory % track_memory(kind(this % ij), capacity * 2, 0, 'MATRIXCACHE::MATRIXARRAY::IJ')
277  call master_memory % track_memory(kind(this % coefficient), capacity, 0, 'MATRIXCACHE::MATRIXARRAY::COEFF')
278  allocate(this % ij(2, capacity), this % coefficient(capacity), stat = err)
279  this % constructed = .true.
280  end if
281 
282  end subroutine construct_array
283 
284 
285  subroutine clear_array (this)
286  class(matrixarray) :: this
287 
288  this % num_elems = 0
289 
290  end subroutine clear_array
291 
292 
293  subroutine insert_into_array (this, i, j, coeff)
294  class(matrixarray) :: this
295  integer, intent(in) :: i, j
296  real(wp), intent(in) :: coeff
297 
298  this % num_elems = this % num_elems + 1
299  this % ij(this % num_elems, 1) = i
300  this % ij(this % num_elems, 2) = j
301  this % coefficient(this % num_elems) = coeff
302 
303  end subroutine insert_into_array
304 
305 
306  subroutine get_from_array(this, idx, i, j, coeff)
307  class(matrixarray), intent(in) :: this
308  integer, intent(in) :: idx
309  integer, intent(out) :: i, j
310  real(wp), intent(out) :: coeff
311 
312  if (idx > this % num_elems) stop "Matrix Array access segfault!"
313 
314  i = this % ij(1, idx)
315  j = this % ij(2, idx)
316  coeff = this % coefficient(idx)
317 
318  end subroutine get_from_array
319 
320 
321  subroutine expand_capacity (this, capacity)
322  class(matrixcache) :: this
323  integer, intent(in) :: capacity
324 
325  do while (this % max_capacity < capacity)
326  call this % expand_array
327  end do
328 
329  end subroutine expand_capacity
330 
331 
332  logical function check_bounds (this, i)
333  class(matrixcache), intent(in) :: this
334  integer, intent(in) :: i
335 
336  if (i <= 0 .or. i > this % n) then
337  write (stdout, "('MatrixCache::check_bounds - Out of Bounds access', 2i4)") i, this % n
338  stop "MatrixCache::check_bounds - Out of Bounds access"
339  check_bounds = .false.
340  else
341  check_bounds = .true.
342 
343  end if
344 
345  end function check_bounds
346 
347 
348  logical function is_empty (this)
349  class(matrixcache), intent(in) :: this
350 
351  is_empty = (this % n == 0)
352 
353  end function is_empty
354 
355 
356  subroutine clear (this)
357  class(matrixcache) :: this
358  integer :: ido
359 
360  if (allocated(this % matrix_arrays)) this % matrix_arrays(:) % num_elems = 0
361  this % n = 0
362 
363  end subroutine clear
364 
365 
366  integer function get_size (this)
367  class(matrixcache), intent(in) :: this
368 
369  get_size = this % n
370 
371  end function get_size
372 
373 
374  subroutine shrink_capacity (this)
375  class(matrixcache) :: this
376  integer :: ido, total_size, shrink_size
377 
378  type(matrixarray), allocatable :: temp_matrix(:)
379 
380  if (allocated(this % matrix_arrays)) then
381  total_size = size(this % matrix_arrays)
382  shrink_size = min(this % n / this % expand_size + 1, this % num_matrix_chunks)
383 
384  if (shrink_size == this % num_matrix_chunks) return
385 
386  do ido = shrink_size + 1, total_size
387  call this % matrix_arrays(ido) % destroy
388  end do
389 
390  allocate(temp_matrix(shrink_size))
391  call master_memory % track_memory(36, size(temp_matrix), 0, 'MATRIXCACHE::SHRINK::TEMP')
392  temp_matrix(1:shrink_size) = this % matrix_arrays(1:shrink_size)
393  call master_memory % free_memory(36, size(this % matrix_arrays))
394  deallocate(this % matrix_arrays)
395 
396  allocate(this % matrix_arrays(shrink_size))
397  call master_memory % track_memory(36, size(this % matrix_arrays), 0, 'MATRIXCACHE::SHRINK::MATRIXARRAYS')
398  this % matrix_arrays(1:shrink_size) = temp_matrix(1:shrink_size)
399  this % max_capacity = shrink_size * this % expand_size
400  this % num_matrix_chunks = shrink_size
401  call master_memory % free_memory(36, size(temp_matrix))
402  deallocate(temp_matrix)
403  end if
404 
405  end subroutine shrink_capacity
406 
407 
408  subroutine clear_and_shrink (this)
409  class(matrixcache) :: this
410 
411  call this % clear
412  call this % shrink_capacity
413 
414  end subroutine clear_and_shrink
415 
416 
417  subroutine destroy_cache (this)
418  class(matrixcache) :: this
419  integer :: num_arrays, ido
420 
421  if (allocated(this % matrix_arrays)) then
422  num_arrays = size(this % matrix_arrays)
423 
424  do ido = 1, num_arrays
425  call this % matrix_arrays(ido) % destroy
426  end do
427 
428  call master_memory % free_memory(36, size(this % matrix_arrays))
429  deallocate(this % matrix_arrays)
430  end if
431 
432  this % constructed = .false.
433 
434  end subroutine destroy_cache
435 
436 
437  subroutine print_matelems (this)
438  class(matrixcache) :: this
439  integer :: labels(8), ido, arrs, jdo, elm
440  real(wp) :: coeff
441 
442  write (stdout, "('Outputting Matrix elements....')")
443 
444  this % matrix_index = this % n / this % expand_size + 1
445  elm = 0
446 
447  do ido = 1, this % num_matrix_chunks
448  write (stdout, *) this % matrix_arrays(ido) % num_elems
449  do jdo = 1, this % matrix_arrays(ido) % num_elems
450  write (stdout, "(2i8,' -- ',D16.8)") this % matrix_arrays(ido) % ij(1:2,jdo), &
451  this % matrix_arrays(ido) % coefficient(jdo)
452  end do
453  end do
454 
455  end subroutine print_matelems
456 
457 
458  subroutine qsort_matelem (this)
459  class(matrixcache) :: this
460 
461  call this % check_constructed
462  if (this % n <= 1) return
463  call this % quick_sort_function(1, this % n)
464  !call QsortMatrixElement_ji(this % matrix_elements(1:this % n))
465 
466  end subroutine qsort_matelem
467 
468 
469  recursive subroutine quick_sort_function (this, start_idx, end_idx)
470  class(matrixcache) :: this
471  integer, intent(in) :: start_idx, end_idx
472  integer :: mat_size, iq
473 
474  mat_size = end_idx - start_idx + 1
475 
476  if (mat_size > 1) then
477  call this % partition_function(iq, start_idx, end_idx)
478  call this % quick_sort_function(start_idx, iq - 1)
479  call this % quick_sort_function(iq, end_idx)
480  end if
481 
482  end subroutine quick_sort_function
483 
484 
485  subroutine partition_function (this, marker, start_idx, end_idx)
486  class(matrixcache) :: this
487  integer, intent(in) :: start_idx, end_idx
488  integer, intent(out) :: marker
489  integer :: i, j, temp_ij(2), x_ij(2), A_ij(2), mat_id_i, mat_id_j, local_index_i, local_index_j
490  real(wp) :: temp_coeff, x_coeff, A_coeff
491 
492  call this % get_from_cache(start_idx, x_ij(1), x_ij(2), x_coeff)
493 
494  i = start_idx - 1
495  j = end_idx + 1
496 
497  do
498  j = j - 1
499  do
500  call this % get_from_cache(j, a_ij(1), a_ij(2), a_coeff)
501  if (a_ij(1) <= x_ij(1)) exit
502  j = j - 1
503  end do
504  i = i + 1
505  do
506  call this % get_from_cache(i, a_ij(1), a_ij(2), a_coeff)
507  if (a_ij(1) >= x_ij(1)) exit
508  i = i + 1
509  end do
510  if (i < j) then
511  ! exchange A(i) and A(j)
512  call this % get_local_mat(i, mat_id_i, local_index_i)
513  call this % get_local_mat(j, mat_id_j, local_index_j)
514 
515  temp_ij(:) = this % matrix_arrays(mat_id_i) % ij(:,local_index_i)
516  temp_coeff = this % matrix_arrays(mat_id_i) % coefficient(local_index_i)
517 
518  this % matrix_arrays(mat_id_i) % ij(:,local_index_i) &
519  = this % matrix_arrays(mat_id_j) % ij(:,local_index_j)
520  this % matrix_arrays(mat_id_i) % coefficient(local_index_i) &
521  = this % matrix_arrays(mat_id_j) % coefficient(local_index_j)
522 
523  this % matrix_arrays(mat_id_j) % ij(:,local_index_j) = temp_ij(:)
524  this % matrix_arrays(mat_id_j) % coefficient(local_index_j) = temp_coeff
525  else if (i == j) then
526  marker = i + 1
527  return
528  else
529  marker = i
530  return
531  end if
532  end do
533 
534  end subroutine partition_function
535 
536  ! Recursive Fortran 95 quicksort routine adapted to
537  ! sort matrix elements into ascending either ij or ji
538  ! Author: Juli Rew, SCD Consulting (juliana@ucar.edu), 9/03
539  ! Based on algorithm from Cormen et al., Introduction to Algorithms,
540  ! 1997 printing
541 
542  ! Made F conformant by Walt Brainerd
543 
544  !subroutine construct(this,i,j,k,l,inttype,positron,coefficent)
545  ! class(MatrixElement),intent(inout) :: this
546  ! integer,intent(in) :: i,j,k,l,positron,inttype
547  ! real(wp),intent(in) :: coefficent
548  !
549  ! this%coefficient = coefficient
550  ! call pack4ints64(i,j,k,l,this%integral_label)
551  ! this%positron = positron
552  !
553  ! if(inttype /= ONE_ELECTRON_INTEGRAL .or. inttype /= TWO_ELECTRON_INTEGRAL) then
554  ! write(stdout,"('Error integral type of ',i4,' is not 1 or 2 electron')") inttype
555  ! stop "Integral label is neither 1 or 2 electron"
556  ! else
557  ! this%integral_type = inttype
558  ! endif
559  !
560  !
561  !
562  !
563  !end subroutine
564 
565 end module matrixcache_module
matrixcache_module::expand_capacity
subroutine expand_capacity(this, capacity)
Definition: MatrixCache_module.f90:322
matrixcache_module::get_chunk_idx
integer function get_chunk_idx(this, idx)
Definition: MatrixCache_module.f90:107
memorymanager_module
Memory manager module.
Definition: MemoryManager_module.f90:31
timing_module
Timing module.
Definition: Timing_Module.f90:28
matrixcache_module::quick_sort_function
recursive subroutine quick_sort_function(this, start_idx, end_idx)
Definition: MatrixCache_module.f90:470
matrixcache_module::destroy_matrix_array
subroutine destroy_matrix_array(this)
Definition: MatrixCache_module.f90:251
matrixcache_module::clear_and_shrink
subroutine clear_and_shrink(this)
Definition: MatrixCache_module.f90:409
matrixcache_module::get_local_mat
subroutine get_local_mat(this, idx, mat_id, local_index)
Definition: MatrixCache_module.f90:185
matrixcache_module::print_matelems
subroutine print_matelems(this)
Definition: MatrixCache_module.f90:438
consts_mpi_ci::default_expand_size
integer, parameter default_expand_size
Default expansion size for the cache.
Definition: consts_mpi_ci.f90:133
matrixcache_module::is_empty
logical function is_empty(this)
Definition: MatrixCache_module.f90:349
matrixcache_module::get_from_cache
subroutine get_from_cache(this, idx, i, j, coeff)
Definition: MatrixCache_module.f90:203
matrixcache_module::get_size
integer function get_size(this)
Definition: MatrixCache_module.f90:367
matrixcache_module::shrink_capacity
subroutine shrink_capacity(this)
Definition: MatrixCache_module.f90:375
matrixcache_module::get_from_array
subroutine get_from_array(this, idx, i, j, coeff)
Definition: MatrixCache_module.f90:307
matrixcache_module::clear
subroutine clear(this)
Definition: MatrixCache_module.f90:357
consts_mpi_ci
MPI SCATCI Constants module.
Definition: consts_mpi_ci.f90:31
matrixcache_module::construct
subroutine construct(this, expand_size_)
Definition: MatrixCache_module.f90:119
matrixcache_module::check_constructed
subroutine check_constructed(this)
Definition: MatrixCache_module.f90:96
matrixcache_module::expand_array
subroutine expand_array(this)
Definition: MatrixCache_module.f90:223
matrixcache_module::matrixarray
Definition: MatrixCache_module.f90:42
matrixcache_module::check_bounds
logical function check_bounds(this, i)
Definition: MatrixCache_module.f90:333
matrixcache_module::partition_function
subroutine partition_function(this, marker, start_idx, end_idx)
Definition: MatrixCache_module.f90:486
matrixcache_module::qsort_matelem
subroutine qsort_matelem(this)
Definition: MatrixCache_module.f90:459
matrixcache_module::insert_into_array
subroutine insert_into_array(this, i, j, coeff)
Definition: MatrixCache_module.f90:294
matrixcache_module::construct_array
subroutine construct_array(this, capacity)
Definition: MatrixCache_module.f90:271
matrixcache_module
Matrix cache module.
Definition: MatrixCache_module.f90:28
matrixcache_module::destroy_cache
subroutine destroy_cache(this)
Definition: MatrixCache_module.f90:418
timing_module::master_timer
type(timer), public master_timer
Definition: Timing_Module.f90:74
matrixcache_module::insert_into_cache
subroutine insert_into_cache(this, i, j, coefficient)
Definition: MatrixCache_module.f90:160
memorymanager_module::master_memory
type(memorymanager), public master_memory
This is the global memory manager.
Definition: MemoryManager_module.f90:69
matrixcache_module::clear_array
subroutine clear_array(this)
Definition: MatrixCache_module.f90:286
matrixcache_module::matrixcache
This handles the matrix elements and also expands the vector size if we have reached max capacity.
Definition: MatrixCache_module.f90:60