MPI-SCATCI  2.0
An MPI version of SCATCI
ELPAMatrix_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 blas_lapack, only: blasint
31  use parallelization_module, only: grid => process_grid
32  use precisn, only: wp
34  use timing_module, only: master_timer
35 
36  implicit none
37 
38  private
39 
47  type, public, extends(scalapackmatrix) :: elpamatrix
48  contains
49  procedure, public :: finalize_matrix => finalize_elpa
50  end type
51 
52 contains
53 
54  subroutine finalize_elpa (this)
55 
56  class(elpamatrix) :: this
57  integer(blasint) :: one = 1, i, k, l, iprow, ipcol
58  real(wp) :: alpha = 1
59 
60  ! finalize parent
61  call this % SCALAPACKMatrix % finalize_matrix
62 
63  ! perform symmetrization (add transposed triangle)
64  call master_timer % start_timer("Symmetrize matrix")
65  call pdtradd('U', 'T', this % mat_dimen, this % mat_dimen, &
66  alpha, this % a_local_matrix, one, one, this % descr_a_mat, &
67  alpha, this % a_local_matrix, one, one, this % descr_a_mat, this % a_local_matrix)
68 
69  ! scale the diagonal
70  do i = 1, this % mat_dimen
71  call infog2l(i, i, this % descr_a_mat, grid % gprows, grid % gpcols, grid % mygrow, grid % mygcol, k, l, iprow, ipcol)
72  if (iprow == grid % mygrow .and. ipcol == grid % mygcol) then
73  this % a_local_matrix(k, l) = this % a_local_matrix(k, l) / 2
74  end if
75  end do
76  call master_timer % stop_timer("Symmetrize matrix")
77 
78  end subroutine finalize_elpa
79 
80 end module elpamatrix_module
timing_module
Timing module.
Definition: Timing_Module.f90:28
elpamatrix_module
ELPA matrix module.
Definition: ELPAMatrix_module.f90:28
parallelization_module::process_grid
type(processgrid) process_grid
Definition: Parallelization_module.F90:80
parallelization_module
Distribution of processes into a grid.
Definition: Parallelization_module.F90:29
scalapackmatrix_module::scalapackmatrix
Definition: SCALAPACKMatrix_module.f90:47
elpamatrix_module::finalize_elpa
subroutine finalize_elpa(this)
Definition: ELPAMatrix_module.f90:55
timing_module::master_timer
type(timer), public master_timer
Definition: Timing_Module.f90:74
scalapackmatrix_module
ScaLAPACK matrix module.
Definition: SCALAPACKMatrix_module.f90:32
elpamatrix_module::elpamatrix
ELPA distributed matrix.
Definition: ELPAMatrix_module.f90:47