MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
DiagonalizerResult_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 Diagonalizer result data object
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
27!>
28module diagonalizerresult_module
29
30 use cdenprop_defs, only: civect
32 use baseintegral_module, only: baseintegral
33 use options_module, only: options
34
35 implicit none
36
37 public diagonalizerresult
38
39 private
40
41 !> \brief Output from diagonalization
42 !> \authors A Al-Refaie
43 !> \date 2017
44 !>
45 !> Base class for \ref SolutionHandler_module::SolutionHandler and \ref Target_RMat_CI_module::Target_RMat_CI.
46 !> Can be used as a storage for the eigensolution when that is to be passed to CDENPROP.
47 !>
48 type, abstract :: diagonalizerresult
49 type(civect) :: ci_vec !< Data structure for keeping the solution of the eigensystem.
50
51 integer :: vector_storage = save_continuum_coeffs
52 integer :: continuum_dimen
53 contains
54 procedure(handle_eigenvalues), deferred, public :: handle_eigenvalues
55 procedure(handle_eigenvector), deferred, public :: handle_eigenvector
56 procedure(export_eigenvalues), deferred, public :: export_eigenvalues
57 procedure :: export_header
58 procedure :: write_header
59 procedure :: finalize_solutions
60 end type diagonalizerresult
61
62 abstract interface
63 !> \brief Main build routine of the hamiltonian
64 !> \authors A Al-Refaie
65 !> \date 2017
66 !>
67 !> All build must be done within this routine in order to be used by MPI-SCATCI.
68 !>
69 !> \param[out] matrix_elements Resulting matrix elements from the build.
70 !>
71 subroutine handle_eigenvalues (this, eigenvalues, diagonals, num_eigenpairs, vec_dimen)
72 use precisn, only : wp
73 import :: diagonalizerresult
74 class(diagonalizerresult) :: this
75 integer, intent(in) :: num_eigenpairs, vec_dimen
76 real(wp), intent(in) :: eigenvalues(num_eigenpairs), diagonals(vec_dimen)
77 end subroutine handle_eigenvalues
78
79 !> \brief Main build routine of the hamiltonian
80 !> \authors A Al-Refaie
81 !> \date 2017
82 !>
83 !> All build must be done within this routine in order to be used by MPI-SCATCI.
84 !>
85 !> \param[out] matrix_elements Resulting matrix elements from the build.
86 !>
87 subroutine handle_eigenvector (this, eigenvector, vec_dimen)
88 use precisn, only : wp
89 import :: diagonalizerresult
90 class(diagonalizerresult) :: this
91 integer, intent(in) :: vec_dimen
92 real(wp), intent(in) :: eigenvector(vec_dimen)
93 end subroutine handle_eigenvector
94
95 !> \brief build routine of the hamiltonian
96 !> \authors A Al-Refaie
97 !> \date 2017
98 !>
99 !> All build must be done within this routine in order to be used by MPI-SCATCI.
100 !>
101 !> \param[out] matrix_elements Resulting matrix elements from the build.
102 !>
103 subroutine export_eigenvalues (this, eigenvalues, diagonals, num_eigenpairs, vec_dimen, ei, iphz)
104 use precisn, only : wp
105 import :: diagonalizerresult
106 class(diagonalizerresult) :: this
107 integer, intent(in) :: num_eigenpairs, vec_dimen
108 real(wp), intent(in) :: eigenvalues(num_eigenpairs), diagonals(vec_dimen)
109 real(wp), allocatable :: ei(:)
110 integer, allocatable :: iphz(:)
111 end subroutine export_eigenvalues
112 end interface
113
114contains
115
116 !> \brief Save eigenvector information internally
117 !> \authors J Benda
118 !> \date 2019
119 !>
120 !> This dummy function is expected to be overridden by derived types that support passing
121 !> the eigensolution to CDENPROP. This is used in \ref SolutionHandler_module::SolutionHandler
122 !> after the main diagonalization.
123 !>
124 subroutine export_header (this, option, integrals)
125 class(diagonalizerresult) :: this
126 class(options), intent(in) :: option
127 class(baseintegral), intent(in) :: integrals
128 end subroutine export_header
129
130
131 !> \brief Save eigenvector information to disk
132 !> \authors J Benda
133 !> \date 2019
134 !>
135 !> This dummy function is expected to be overridden by derived types that support writing
136 !> the eigensolution to disk file. This is used in \ref SolutionHandler_module::SolutionHandler
137 !> after the main diagonalization.
138 !>
139 subroutine write_header (this, option, integrals)
140 class(diagonalizerresult) :: this
141 class(options), intent(in) :: option
142 class(baseintegral), intent(in) :: integrals
143 end subroutine write_header
144
145
146 !> \brief Finalize saving eigenvector information to disk
147 !> \authors J Benda
148 !> \date 2019
149 !>
150 !> This dummy function is expected to be overridden by derived types that support writing
151 !> the eigensolution to disk file. This is used in \ref SolutionHandler_module::SolutionHandler
152 !> after the main diagonalization.
153 !>
154 subroutine finalize_solutions (this, option)
155 class(diagonalizerresult) :: this
156 class(options), intent(in) :: option
157 end subroutine finalize_solutions
158
159end module diagonalizerresult_module
MPI SCATCI Constants module.
integer, parameter save_continuum_coeffs
Omit L2 section of the eigenvectors (only keep continuum channels)