MPI-SCATCI  2.0
An MPI version of SCATCI
ARPACKDiagonalizer_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 precisn, only: longint, wp
34  use const_gbl, only: stdout
35  use basematrix_module, only: basematrix
40  use options_module, only: options
41  use mpi_gbl, only: myrank, master, mpi_xermsg
42  use m_4_mkarp, only: mkarp
44 
45  implicit none
46 
48  contains
49  procedure, public :: diagonalize => diagonalize_arpack
50  procedure, private :: diagonalize_writermatrix
51  end type arpackdiagonalizer
52 
53 contains
54 
55  subroutine diagonalize_arpack (this, matrix_elements, num_eigenpair, dresult, all_procs, option, integrals)
56  class(arpackdiagonalizer) :: this
57  class(diagonalizerresult) :: dresult
58  class(basematrix), intent(in) :: matrix_elements
59  class(baseintegral), intent(in) :: integrals
60  type(options), intent(in) :: option
61  integer, intent(in) :: num_eigenpair
62  logical, intent(in) :: all_procs
63 
64  integer :: max_iterations
65  real(wp) :: max_tolerance
66 
67  !Only the first rank can do it
68  if (myrank /= master) return
69 
70  if (iand(dresult % vector_storage, pass_to_cdenprop) /= 0) then
71  call mpi_xermsg('ARPACKDiagonalizer_module', 'diagonalize_ARPACK', &
72  'PASS_TO_CDENPROP not implemented for Arpack', 1, 1)
73  end if
74 
75  max_iterations = option % max_iterations
76  max_tolerance = option % max_tolerance
77 
78  if (max_iterations < 0) then
79  max_iterations = max(num_eigenpair * 50, 500)
80  end if
81 
82  select type(matrix_elements)
83  type is (writermatrix)
84  call this % diagonalize_writermatrix(matrix_elements, num_eigenpair, dresult, max_iterations, max_tolerance, &
85  option, integrals)
86  class is (basematrix)
87  call mpi_xermsg('ARPACKDiagonalizer_module', 'diagonalize_ARPACK', &
88  'Only WriterMatrix format can use Arpack', 1, 1)
89  end select
90 
91  end subroutine diagonalize_arpack
92 
93 
94  subroutine diagonalize_writermatrix (this, matrix_elements, num_eigenpair, dresult, max_iterations, max_tolerance, &
95  option, integrals)
96  class(arpackdiagonalizer) :: this
97  class(diagonalizerresult) :: dresult
98  type(writermatrix), intent(in) :: matrix_elements
99  class(baseintegral), intent(in) :: integrals
100  type(options), intent(in) :: option
101  integer, intent(in) :: num_eigenpair
102  integer, intent(in) :: max_iterations
103  real(wp), intent(in) :: max_tolerance
104 
105  !for now these are default but will be included i the input
106 
107  integer :: matrix_unit
108  integer :: matrix_size
109  integer :: i1, i2
110 
111  integer :: arpack_numeig
112  integer :: arpack_maxit
113  real(wp) :: arpack_maxtol
114  integer :: num_matrix_elements_per_record, num_elems
115 
116  real(wp), allocatable :: eigenvalues(:), eigenvector(:)
117 
118  write (stdout, "('Diagonalization done with ARPACK')")
119  write (stdout, "('Parameters:')")
120  write (stdout, "('N: ',i8)") matrix_elements % get_matrix_size()
121  write (stdout, "('Requested # of eigenpairs',i8)") num_eigenpair
122 
123  arpack_numeig = num_eigenpair
124  arpack_maxit = max_iterations
125  arpack_maxtol = max_tolerance
126  matrix_size = matrix_elements % get_matrix_size()
127  num_elems = matrix_elements % get_size()
128  matrix_unit = matrix_elements % get_matrix_unit()
129 
130  allocate(eigenvalues(num_eigenpair), eigenvector(matrix_size))
131 
132  call mkarp(matrix_size, arpack_numeig, arpack_maxit, arpack_maxtol, matrix_unit)
133 
134  ! 3. Opening the produced file "eigensystem_arpack_file.bin" file.
135  ! This file contains the produced eigenpairs and it is
136  ! written to the hard disk.
137 
138  open (unit = 11, file = "___tmp_eigensys", status = "old", action = "read", form = "unformatted")
139 
140  ! 4. Reading the eigenvalues.
141 
142  do i1 = 1, num_eigenpair
143  read (unit = 11) eigenvalues(i1)
144  end do
145 
146  if (iand(option % vector_storage_method, pass_to_cdenprop) /= 0) then
147  call dresult % export_header(option, integrals)
148  end if
149  call dresult % write_header(option, integrals)
150  call dresult % handle_eigenvalues(eigenvalues, matrix_elements % diagonal, num_eigenpair, matrix_size)
151 
152  ! 5. Reading the eigenvectors.
153 
154  do i1 = 1, num_eigenpair
155  do i2 = 1, matrix_size
156  read (unit = 11) eigenvector(i2)
157  end do
158  call dresult % handle_eigenvector(eigenvector, matrix_size)
159  eigenvector = 0
160  end do
161 
162  ! 6. Closing the unit and deleting the file
163  ! "___tmp_eigensys".
164 
165  close (unit = 11, status = "delete")
166 
167  deallocate(eigenvalues, eigenvector)
168 
169  end subroutine diagonalize_writermatrix
170 
171 end module arpackdiagonalizer_module
writermatrix_module
Write matrix module.
Definition: WriterMatrix_module.f90:31
diagonalizerresult_module::diagonalizerresult
Output from diagonalization.
Definition: DiagonalizerResult_module.f90:48
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
arpackdiagonalizer_module
Diagonalizer type using Arpack backend.
Definition: ARPACKDiagonalizer_module.f90:31
arpackdiagonalizer_module::diagonalize_arpack
subroutine diagonalize_arpack(this, matrix_elements, num_eigenpair, dresult, all_procs, option, integrals)
Definition: ARPACKDiagonalizer_module.f90:56
basematrix_module
Base matrix module.
Definition: BaseMatrix_module.f90:28
consts_mpi_ci
MPI SCATCI Constants module.
Definition: consts_mpi_ci.f90:31
arpackdiagonalizer_module::arpackdiagonalizer
Definition: ARPACKDiagonalizer_module.f90:47
basematrix_module::basematrix
Base matrix type.
Definition: BaseMatrix_module.f90:47
baseintegral_module
Base integral module.
Definition: BaseIntegral_module.f90:28
diagonalizerresult_module
Diagonalizer result data object.
Definition: DiagonalizerResult_module.f90:28
diagonalizer_module
Base type for all diagonalizers.
Definition: Diagonalizer_module.f90:30
diagonalizer_module::basediagonalizer
Definition: Diagonalizer_module.f90:34
baseintegral_module::baseintegral
Definition: BaseIntegral_module.f90:41
arpackdiagonalizer_module::diagonalize_writermatrix
subroutine diagonalize_writermatrix(this, matrix_elements, num_eigenpair, dresult, max_iterations, max_tolerance, option, integrals)
Definition: ARPACKDiagonalizer_module.f90:96
options_module
Options module.
Definition: Options_module.f90:41
consts_mpi_ci::pass_to_cdenprop
integer, parameter pass_to_cdenprop
Eigenpairs will be saved in memory and passed to CDENPROP and outer interface.
Definition: consts_mpi_ci.f90:110