MPI-SCATCI  2.0
An MPI version of SCATCI
WriterMatrix_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 
32 
33  use const_gbl, only: stdout
34  use mpi_gbl, only: master
35  use precisn, only: wp
36  use scatci_routines, only: wrtem
39  use options_module, only: options
40  use parallelization_module, only: grid => process_grid
41 
42  implicit none
43 
44  public writermatrix
45 
46  private
47 
56  type, extends(distributedmatrix) :: writermatrix
59  type(matrixcache) :: write_cache
60 
61  integer :: num_elements_per_record
62  integer :: matrix_unit
63  contains
64  procedure, public :: print => print_writer
65  procedure, public :: finalize_matrix_self => finalize_matrix_writer
66  procedure, public :: set_options => set_options_writer
67  procedure, public :: get_matelem_self => get_matelem_writer
68  procedure, public :: clear_matrix => clear_writer
69  procedure, public :: destroy_matrix => destroy_writer
70  procedure, public :: write_cache_to_unit
71  procedure, public :: insert_into_diag_matrix => insert_into_write_cache
72  procedure, public :: get_matrix_unit
73  end type writermatrix
74 
75 contains
76 
81  integer function get_matrix_unit (this)
82  class(writermatrix) :: this
83 
84  get_matrix_unit = this % matrix_unit
85 
86  end function get_matrix_unit
87 
88 
93  subroutine set_options_writer (this, option_val)
94  class(writermatrix) :: this
95  class(options), intent(in) :: option_val
96 
97  this % matrix_unit = option_val % hamiltonian_unit
98  this % num_elements_per_record = option_val % num_matrix_elements_per_rec
99 
100  !we set the expansion size to the number of records + 1
101  call this % write_cache % construct(this % num_elements_per_record + 1)
102 
103  end subroutine set_options_writer
104 
105 
110  subroutine get_matelem_writer (this, idx, i, j, coeff)
111  class(writermatrix) :: this
112  integer, intent(in) :: idx
113  integer, intent(out) :: i, j
114  real(wp), intent(out) :: coeff
115 
116  i = -1
117  j = -1
118  coeff = -1.0
119 
120  stop "get element not yet implemented"
121 
122  end subroutine get_matelem_writer
123 
124 
132  subroutine write_cache_to_unit (this)
133  class(writermatrix) :: this
134  integer :: record_count
135 
136  if (this % write_cache % is_empty()) return
137 
138  record_count = this % write_cache % get_size()
139 
140  !Write it all out
141  call wrtem(this % matrix_unit, record_count, this % num_elements_per_record, &
142  this % write_cache % matrix_arrays(1) % ij(1:2,1:record_count), &
143  this % write_cache % matrix_arrays(1) % coefficient(1:record_count))
144  call this % write_cache % clear
145 
146  end subroutine write_cache_to_unit
147 
148 
158  logical function insert_into_write_cache (this, row, column, coefficient)
159  class(writermatrix) :: this
160  integer, intent(in) :: row, column
161  real(wp), intent(in) :: coefficient
162  integer :: record_count
163 
164  if (grid % grank /= master) then
165  insert_into_write_cache = .false.
166  return
167  end if
168 
169  if (row == column) then
170  if (row > size(this % diagonal)) then
171  write (stdout, "('Row greater than size of diagonal!!!!',2i12)") row, size(this % diagonal)
172  stop "Error"
173  end if
174  this % diagonal(row) = coefficient
175  end if
176 
177  insert_into_write_cache = .true.
178  if (abs(coefficient) < this % threshold) return
179 
180  call this % write_cache % insert_into_cache(row, column, coefficient)
181  record_count = this % write_cache % get_size()
182 
183  !Update number of elements
184  this % n = this % n + 1
185  insert_into_write_cache = .true.
186 
187  if (record_count == this % num_elements_per_record) then
188  !Write it all out
189  call this % write_cache_to_unit
190  end if
191 
192  end function insert_into_write_cache
193 
194 
202  subroutine finalize_matrix_writer (this)
203  class(writermatrix) :: this
204  integer :: dummy_ij(2, this % num_elements_per_record), record_count = 0
205  real(wp) :: dummy_coeff(this % num_elements_per_record)
206 
207  if (grid % grank /= master) return
208 
209  write (stdout, "('finishing unit')")
210  call this % write_cache_to_unit
211  record_count = 0
212  call wrtem(this % matrix_unit, record_count, this % num_elements_per_record, dummy_ij, dummy_coeff)
213  write (stdout, "('done')")
214 
215  end subroutine finalize_matrix_writer
216 
217 
222  subroutine print_writer (this)
223  class(writermatrix) :: this
224 
225  write (stdout, "('-------TEMP CACHE---------')")
226  call this % temp_cache % print
227  !write(stdout,"('-------HARD CACHE---------')")
228  !call this%matrix_cache%print
229 
230  end subroutine print_writer
231 
232 
239  subroutine clear_writer (this)
240  class(writermatrix) :: this
241 
242  call this % write_cache % clear
243 
244  end subroutine clear_writer
245 
246 
251  subroutine destroy_writer (this)
252  class(writermatrix) :: this
253 
254  call this % clear
255  call this % write_cache % destroy
256 
257  end subroutine destroy_writer
258 
259 end module writermatrix_module
writermatrix_module
Write matrix module.
Definition: WriterMatrix_module.f90:31
options_module::options
Calculation setup for a single Hamiltonian diagonalization.
Definition: Options_module.f90:74
writermatrix_module::writermatrix
Matrix associated with a disk drive.
Definition: WriterMatrix_module.f90:56
writermatrix_module::clear_writer
subroutine clear_writer(this)
Clear matrix.
Definition: WriterMatrix_module.f90:240
writermatrix_module::insert_into_write_cache
logical function insert_into_write_cache(this, row, column, coefficient)
Insert matrix element.
Definition: WriterMatrix_module.f90:159
writermatrix_module::destroy_writer
subroutine destroy_writer(this)
Destroy matrix.
Definition: WriterMatrix_module.f90:252
distributedmatrix_module::distributedmatrix
Distributed matrix type.
Definition: DistributedMatrix_module.f90:57
writermatrix_module::write_cache_to_unit
subroutine write_cache_to_unit(this)
Write matrix chunk to file.
Definition: WriterMatrix_module.f90:133
parallelization_module::process_grid
type(processgrid) process_grid
Definition: Parallelization_module.F90:80
distributedmatrix_module
Distributed matrix module.
Definition: DistributedMatrix_module.f90:32
parallelization_module
Distribution of processes into a grid.
Definition: Parallelization_module.F90:29
writermatrix_module::print_writer
subroutine print_writer(this)
Print matrix to standard output.
Definition: WriterMatrix_module.f90:223
writermatrix_module::set_options_writer
subroutine set_options_writer(this, option_val)
Initialize the matrix cache using the Options object.
Definition: WriterMatrix_module.f90:94
writermatrix_module::finalize_matrix_writer
subroutine finalize_matrix_writer(this)
Write Hamiltonian to disk.
Definition: WriterMatrix_module.f90:203
writermatrix_module::get_matelem_writer
subroutine get_matelem_writer(this, idx, i, j, coeff)
Not implemented.
Definition: WriterMatrix_module.f90:111
matrixcache_module
Matrix cache module.
Definition: MatrixCache_module.f90:28
writermatrix_module::get_matrix_unit
integer function get_matrix_unit(this)
Return Hamiltonian disk file unit.
Definition: WriterMatrix_module.f90:82
options_module
Options module.
Definition: Options_module.f90:41
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