MPI-SCATCI  2.0
An MPI version of SCATCI
ELPADiagonalizer_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 
32  use blas_lapack, only: blasint
33  use const, only: stdout
34  use mpi_mod, only: mpi_xermsg, nprocs, myrank
35  use options_module, only: options
36  use parallelization_module, only: grid => process_grid
37  use precisn, only: wp
38  use elpamatrix_module, only: elpamatrix
41 
42  use iso_c_binding, only: c_int
43  use elpa, only: elpa_t, elpa_allocate, elpa_deallocate, elpa_init, elpa_uninit, elpa_ok, elpa_solver_1stage
44 
45  implicit none
46 
47  integer, parameter :: elpaint = c_int
48 
58  contains
59  procedure, public :: diagonalize_backend => diagonalize_backend_elpa
60  end type elpadiagonalizer
61 
62 contains
63 
72  subroutine initialize_elpa
73 
74  integer(elpaint), parameter :: elpa_api = 20181112
75  integer(elpaint) :: success
76 
77  success = elpa_init(elpa_api)
78  call elpa_assert('initialize_elpa', 'ELPA API version not supported', success)
79 
80  end subroutine initialize_elpa
81 
82 
92  subroutine diagonalize_backend_elpa (this, matrix_elements, num_eigenpair, z_mat, w, all_procs, option, integrals)
93 
94  class(elpadiagonalizer) :: this
95  class(scalapackmatrix), intent(in) :: matrix_elements
96  class(baseintegral), intent(in) :: integrals
97  type(options), intent(in) :: option
98  integer, intent(in) :: num_eigenpair
99  logical, intent(in) :: all_procs
100  real(wp), allocatable, intent(inout) :: z_mat(:,:), w(:)
101 
102  integer :: ierr
103  integer(elpaint) :: success, nb, loc_r, loc_c, mat_dimen, myrow, mycol, nev, comm
104  class(elpa_t), pointer :: elpa
105 
106  write (stdout, "('Diagonalization will be done with ELPA')")
107 
108  comm = grid % gcomm
109  myrow = grid % mygrow
110  mycol = grid % mygcol
111  nev = num_eigenpair
112 
113  mat_dimen = matrix_elements % get_matrix_size()
114  loc_r = matrix_elements % local_row_dimen
115  loc_c = matrix_elements % local_col_dimen
116  nb = matrix_elements % scal_block_size
117 
118  ! allocate the eigenvectors
119  allocate(z_mat(loc_r, loc_c), w(mat_dimen), stat = ierr)
120  if (ierr /= 0) then
121  call mpi_xermsg('ELPADiagonalizer_module', 'diagonalize_backend_ELPA', 'Memory allocation error.', ierr, 1)
122  end if
123 
124  ! create the ELPA object
125  elpa => elpa_allocate(success)
126  call elpa_assert('diagonalize_backend_ELPA', 'Failed to create ELPA object.', success)
127 
128  ! set matrix-related ELPA parameters
129  call elpa % set('na', mat_dimen, success)
130  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "na" for ELPA.', success)
131  call elpa % set('nev', mat_dimen, success)
132  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "nev" for ELPA.', success)
133  call elpa % set('local_nrows', loc_r, success)
134  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "local_nrows" for ELPA.', success)
135  call elpa % set('local_ncols', loc_c, success)
136  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "local_ncols" for ELPA.', success)
137  call elpa % set('nblk', nb, success)
138  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "nblk" for ELPA.', success)
139  call elpa % set('mpi_comm_parent', comm, success)
140  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "mpi_comm_parent" for ELPA.', success)
141  call elpa % set('process_row', myrow, success)
142  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "process_row" for ELPA.', success)
143  call elpa % set('process_col', mycol, success)
144  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set "process_col" for ELPA.', success)
145 
146  ! finalize the setup
147  success = elpa % setup()
148  call elpa_assert('diagonalize_backend_ELPA', 'Failed to set up the ELPA object.', success)
149 
150  ! run the solver
151  call elpa % eigenvectors(matrix_elements % a_local_matrix, w, z_mat, success)
152  call elpa_assert('diagonalize_backend_ELPA', 'Diagonalization using ELPA failed.', success)
153 
154  ! free the memory, possible getting ready for another diagonalization
155  call elpa_deallocate(elpa, success)
156  call elpa_assert('diagonalize_backend_ELPA', 'Failed to deallocate the ELPA object.', success)
157 
158  end subroutine diagonalize_backend_elpa
159 
160 
168  subroutine elpa_assert (sub, msg, err)
169 
170  character(len=*), intent(in) :: sub, msg
171  integer(elpaint), intent(in) :: err
172 
173  if (err /= elpa_ok) then
174  call mpi_xermsg('ELPADiagonalizer_module', sub, msg, int(err), 1)
175  end if
176 
177  end subroutine elpa_assert
178 
179 end module elpadiagonalizer_module
elpamatrix_module
ELPA matrix module.
Definition: ELPAMatrix_module.f90:28
options_module::options
Calculation setup for a single Hamiltonian diagonalization.
Definition: Options_module.f90:74
parallelization_module::process_grid
type(processgrid) process_grid
Definition: Parallelization_module.F90:80
baseintegral_module
Base integral module.
Definition: BaseIntegral_module.f90:28
parallelization_module
Distribution of processes into a grid.
Definition: Parallelization_module.F90:29
elpadiagonalizer_module::elpadiagonalizer
Diagonalizer class.
Definition: ELPADiagonalizer_module.f90:57
scalapackmatrix_module::scalapackmatrix
Definition: SCALAPACKMatrix_module.f90:47
elpadiagonalizer_module
Diagonalizer type using ELPA backend.
Definition: ELPADiagonalizer_module.f90:28
diagonalizer_module
Base type for all diagonalizers.
Definition: Diagonalizer_module.f90:30
diagonalizer_module::basediagonalizer
Definition: Diagonalizer_module.f90:34
elpadiagonalizer_module::diagonalize_backend_elpa
subroutine diagonalize_backend_elpa(this, matrix_elements, num_eigenpair, z_mat, w, all_procs, option, integrals)
Diagonalization kernel.
Definition: ELPADiagonalizer_module.f90:93
baseintegral_module::baseintegral
Definition: BaseIntegral_module.f90:41
scalapackdiagonalizer_module::scalapackdiagonalizer
Definition: SCALAPACKDiagonalizer_module.f90:54
elpadiagonalizer_module::elpaint
integer, parameter elpaint
Definition: ELPADiagonalizer_module.f90:47
options_module
Options module.
Definition: Options_module.f90:41
scalapackmatrix_module
ScaLAPACK matrix module.
Definition: SCALAPACKMatrix_module.f90:32
elpamatrix_module::elpamatrix
ELPA distributed matrix.
Definition: ELPAMatrix_module.f90:47
scalapackdiagonalizer_module
Diagonalizer type using ScaLAPACK backend.
Definition: SCALAPACKDiagonalizer_module.f90:32
elpadiagonalizer_module::initialize_elpa
subroutine initialize_elpa
Initialize the ELPA library.
Definition: ELPADiagonalizer_module.f90:73
elpadiagonalizer_module::elpa_assert
subroutine elpa_assert(sub, msg, err)
Check ELPA error code.
Definition: ELPADiagonalizer_module.f90:169