TABLE OF CONTENTS


ABINIT/m_elpa [ Modules ]

[ Top ] [ Modules ]

NAME

 m_elpa

FUNCTION

 This module contains interfaces to ELPA library methods
 See http://elpa.mpcdf.mpg.de

COPYRIGHT

 Copyright (C) 2016-2024 ABINIT group (MT)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

 17 #if defined HAVE_CONFIG_H
 18 #include "config.h"
 19 #endif
 20 
 21 #include "abi_common.h"
 22 
 23 !!#if (defined HAVE_LINALG_ELPA) && (defined HAVE_LINALG_ELPA_2017)
 24 !!#define HAVE_LINALG_ELPA_FORTRAN2008
 25 !!#else
 26 !!#undef HAVE_LINALG_ELPA_FORTRAN2008
 27 !!#endif
 28 
 29 module m_elpa
 30 
 31  use defs_basis
 32  use m_errors
 33 
 34 #ifdef HAVE_LINALG_ELPA
 35 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
 36  use elpa
 37 #else
 38  use elpa1
 39 #endif
 40 #endif
 41 
 42  implicit none
 43 
 44  private
 45 
 46 #ifdef HAVE_LINALG_ELPA
 47 
 48 !Public procedures
 49 !Had to choose names different from those provided by elpa
 50  public :: elpa_func_init                ! Init ELPA
 51  public :: elpa_func_uninit              ! End ELPA
 52  public :: elpa_func_allocate            ! Allocate a ELPA handle and set up MPI information
 53  public :: elpa_func_deallocate          ! Deallocate a ELPA handle
 54  public :: elpa_func_error_handler       ! Manage errors (print a readable message)
 55  public :: elpa_func_get_communicators   ! Get rows and cols communicators (not supposed to be called directly)
 56  public :: elpa_func_set_matrix          ! Set matrix specifications in a ELPA handle
 57  public :: elpa_func_solve_evp_1stage    ! Solve the diagonalization problem (use a ELPA handle)
 58  public :: elpa_func_solve_gevp_2stage   ! Solve the generalized diagonalization problem (use a ELPA handle)
 59  public :: elpa_func_cholesky            ! Apply Cholesky transformation (use a ELPA handle)
 60  public :: elpa_func_invert_triangular   ! Invert triangular matrix (use a ELPA handle)
 61  public :: elpa_func_hermitian_multiply  ! Perform C := A**H * B (use a ELPA handle)
 62 
 63  interface elpa_func_solve_evp_1stage
 64    module procedure elpa_func_solve_evp_1stage_real
 65    module procedure elpa_func_solve_evp_1stage_complex
 66  end interface elpa_func_solve_evp_1stage
 67 
 68  interface elpa_func_solve_gevp_2stage
 69    module procedure elpa_func_solve_gevp_2stage_real
 70    module procedure elpa_func_solve_gevp_2stage_complex
 71  end interface elpa_func_solve_gevp_2stage
 72 
 73  interface elpa_func_cholesky
 74    module procedure elpa_func_cholesky_real
 75    module procedure elpa_func_cholesky_complex
 76  end interface elpa_func_cholesky
 77 
 78  interface elpa_func_invert_triangular
 79    module procedure elpa_func_invert_triangular_real
 80    module procedure elpa_func_invert_triangular_complex
 81  end interface elpa_func_invert_triangular
 82 
 83  interface elpa_func_hermitian_multiply
 84    module procedure elpa_func_hermitian_multiply_real
 85    module procedure elpa_func_hermitian_multiply_complex
 86  end interface elpa_func_hermitian_multiply
 87 
 88 !ELPA generalized handle
 89  type,public :: elpa_hdl_t
 90    logical :: is_allocated=.false.
 91    logical :: matrix_is_set=.false.
 92 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
 93    class(elpa_t),pointer :: elpa
 94 #else
 95    integer :: mpi_comm_parent
 96    integer :: elpa_comm_rows,elpa_comm_cols
 97    integer :: process_row,process_col
 98    integer :: local_nrows,local_ncols
 99    integer :: na,nblk,nev
100    integer :: gpu=0
101 #endif
102  end type elpa_hdl_t
103 
104 #endif
105 
106 CONTAINS  !==============================================================================

m_elpa/elpa_func_allocate [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_allocate

FUNCTION

  Allocate a ELPA handle and set it up with communicators specification

INPUTS

  [blacs_ctx]= -- optional -- Blacs context
  [gpu]= -- optional -- Flag (0 or 1): use GPU version (currently only NVidia)

SIDE EFFECTS

  elpa_hdl(type<elpa_hdl_t>)= ELPA handle

SOURCE

209 subroutine elpa_func_allocate(elpa_hdl,gpu,blacs_ctx)
210 
211 !Arguments ------------------------------------
212  integer,intent(in),optional :: gpu,blacs_ctx
213  type(elpa_hdl_t),intent(inout) :: elpa_hdl
214 
215 !Local variables-------------------------------
216  integer :: err,l_gpu,l_blacs_ctx
217  logical :: gpu_debug_mode=.false.
218  character(len=10) :: varname
219 
220 ! *********************************************************************
221 
222  err=0
223  ! if optional parameter is present, use it
224  ! else use default value, i.e. don't use GPU
225  l_gpu = 0
226  if (present(gpu)) then
227    if(gpu/=0) l_gpu = 1
228  end if
229 
230 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
231  elpa_hdl%elpa => elpa_allocate(err)
232  call elpa_func_error_handler(err_code=err,err_msg='Error in initialization')
233 
234  if(l_gpu==1) then
235 #if defined HAVE_GPU_CUDA
236    varname="nvidia-gpu"
237 #else
238    ABI_BUG("Requested unsupported GPU model with ELPA!")
239 #endif
240    if (err==ELPA_OK) call elpa_hdl%elpa%set(varname,l_gpu,err)
241 
242    ! Handling GPU-support on older ELPA versions (2020.11 and previous)
243    ! Only NVIDIA CUDA GPUs were supported with 'gpu' setting
244    if (err==ELPA_ERROR_ENTRY_NOT_FOUND) then
245 #if defined HAVE_GPU_CUDA
246      varname="gpu"
247      call elpa_hdl%elpa%set(varname,l_gpu,err)
248 #else
249      ABI_ERROR("You seem to use an old version of ELPA ( < 2021.x ) which only supports NVIDIA GPUs.")
250 #endif
251    end if
252    call elpa_func_error_handler(err_code=err,err_msg='Error when enabling GPU on ELPA')
253 
254    if (gpu_debug_mode) then
255      if (err==ELPA_OK) call elpa_hdl%elpa%set("debug",1,err) 
256      call elpa_func_error_handler(err_code=err,err_msg='Error when enabling debug on ELPA')
257    end if
258 
259  end if
260 #else
261  if (err==0.and.l_gpu==1) then
262    elpa_hdl%gpu=l_gpu
263    if (gpu_debug_mode) elpa_hdl%debug=1
264  end if
265 #endif
266 
267  if (present(blacs_ctx)) then
268    if (err==ELPA_OK) call elpa_hdl%elpa%set("blacs_context",int(blacs_ctx,kind=c_int),err)
269    call elpa_func_error_handler(err_code=err,err_varname=varname)
270  end if
271 
272  elpa_hdl%is_allocated=.true.
273 
274 end subroutine elpa_func_allocate

m_elpa/elpa_func_cholesky_complex [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_cholesky_complex

FUNCTION

  Wrapper to elpa_cholesky_complex ELPA function

INPUTS

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix which should be factorized.
                     Distribution is like in Scalapack.
                     Only upper triangle is needs to be set.
                     On return, the upper triangle contains the Cholesky factor
                     and the lower triangle is set to 0.
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

 946 subroutine elpa_func_cholesky_complex(elpa_hdl,aa)
 947 
 948 !Arguments ------------------------------------
 949 !scalars
 950  type(elpa_hdl_t),intent(inout) :: elpa_hdl
 951 !arrays
 952  complex(dpc),intent(inout) :: aa(:,:)
 953 
 954 !Local variables-------------------------------
 955  integer :: err
 956  logical :: success
 957 
 958 ! *********************************************************************
 959 
 960  success=.true. ; err=0
 961 
 962  if (.not.elpa_hdl%is_allocated) then
 963    ABI_BUG('ELPA handle not allocated!')
 964  end if
 965  if (.not.elpa_hdl%matrix_is_set) then
 966    ABI_BUG('Matrix not set in ELPA handle!')
 967  end if
 968 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
 969  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
 970 #else
 971  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
 972 #endif
 973 
 974 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
 975  call elpa_hdl%elpa%cholesky(aa,err)
 976  success=(err==ELPA_OK)
 977 #elif (defined HAVE_LINALG_ELPA_2016)
 978  success = elpa_cholesky_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
 979 &                                elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
 980 #elif (defined HAVE_LINALG_ELPA_2015_11)
 981  call cholesky_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
 982 &                      elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
 983 #elif (defined HAVE_LINALG_ELPA_2015_02)
 984  call cholesky_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
 985 &                      elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
 986 #elif (defined HAVE_LINALG_ELPA_2014)
 987  call cholesky_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
 988 &                      elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,success)
 989 #elif (defined HAVE_LINALG_ELPA_2013)
 990  call cholesky_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
 991 &                      elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
 992 #else
 993 !ELPA-LEGACY-2017
 994  success = elpa_cholesky_complex_double(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
 995 &                                       elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
 996 #endif
 997 
 998  if (.not.success) call elpa_func_error_handler(err_msg='Error in solve_cholesky_complex!')
 999 
1000 end subroutine elpa_func_cholesky_complex

m_elpa/elpa_func_cholesky_real [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_cholesky_real

FUNCTION

  Wrapper to elpa_cholesky_real ELPA function

INPUTS

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix which should be factorized.
                     Distribution is like in Scalapack.
                     Only upper triangle is needs to be set.
                     On return, the upper triangle contains the Cholesky factor
                     and the lower triangle is set to 0.
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

869 subroutine elpa_func_cholesky_real(elpa_hdl,aa)
870 
871 !Arguments ------------------------------------
872 !scalars
873  type(elpa_hdl_t),intent(inout) :: elpa_hdl
874 !arrays
875  real(dp),intent(inout) :: aa(:,:)
876 
877 !Local variables-------------------------------
878  integer :: err
879  logical :: success
880 
881 ! *********************************************************************
882 
883  success=.true. ; err=0
884 
885  if (.not.elpa_hdl%is_allocated) then
886    ABI_BUG('ELPA handle not allocated!')
887  end if
888  if (.not.elpa_hdl%matrix_is_set) then
889    ABI_BUG('Matrix not set in ELPA handle!')
890  end if
891 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
892  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
893 #else
894  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
895 #endif
896 
897 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
898  call elpa_hdl%elpa%cholesky(aa,err)
899  success=(err==ELPA_OK)
900 #elif (defined HAVE_LINALG_ELPA_2016)
901  success = elpa_cholesky_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
902 &                             elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
903 #elif (defined HAVE_LINALG_ELPA_2015_11)
904  call cholesky_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
905                     elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
906 #elif (defined HAVE_LINALG_ELPA_2015_02)
907  call cholesky_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
908                     elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
909 #elif (defined HAVE_LINALG_ELPA_2014)
910  call cholesky_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
911 &                   elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,success)
912 #elif (defined HAVE_LINALG_ELPA_2013)
913  call cholesky_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
914 &                   elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
915 #else
916 !ELPA-LEGACY-2017
917  success = elpa_cholesky_real_double(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
918 &                                    elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
919 #endif
920 
921  if (.not.success) call elpa_func_error_handler(err_msg='Error in solve_cholesky_real!')
922 
923 end subroutine elpa_func_cholesky_real

m_elpa/elpa_func_deallocate [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_deallocate_matrix

FUNCTION

  Deallocate a ELPA handle

INPUTS

SIDE EFFECTS

  elpa_hdl(type<elpa_hdl_t>)= ELPA handle

SOURCE

293 subroutine elpa_func_deallocate(elpa_hdl)
294 
295 !Arguments ------------------------------------
296  type(elpa_hdl_t),intent(inout) :: elpa_hdl
297 
298 !Local variables-------------------------------
299  integer :: err
300 
301 ! *********************************************************************
302 
303  err=0
304 
305 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
306  call elpa_deallocate(elpa_hdl%elpa)
307 #endif
308 
309  elpa_hdl%matrix_is_set=.false.
310  elpa_hdl%is_allocated=.false.
311 
312 end subroutine elpa_func_deallocate

m_elpa/elpa_func_error_handler [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_error_handler

FUNCTION

  Handle ELPA errors

INPUTS

  [err_code]= --optional-- Error code
  [err_msg]= --optional-- Generic error message
  [err_varname]= -- optional-- Name of the ELPA variable related to the error

OUTPUT

  No output, only printing

SOURCE

334 subroutine elpa_func_error_handler(err_code,err_msg,err_varname)
335 
336 !Arguments ------------------------------------
337  integer,optional :: err_code
338  character(len=*),optional :: err_msg,err_varname
339 
340 !Local variables-------------------------------
341  integer :: err_code_
342  character(len=500) :: msg
343  character(len=100) :: err_strg
344 
345 ! *********************************************************************
346 
347  err_code_=-100;if (present(err_code)) err_code_=err_code
348  if (err_code_==0) return
349 
350  err_strg=''
351 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
352  if (err_code_==ELPA_ERROR) err_strg='ELPA_ERROR'
353  if (err_code_==ELPA_ERROR_ENTRY_READONLY) err_strg='ELPA_ERROR_ENTRY_READONLY'
354  if (err_code_==ELPA_ERROR_ENTRY_NOT_FOUND) err_strg='ELPA_ERROR_ENTRY_NOT_FOUND'
355  if (err_code_==ELPA_ERROR_ENTRY_ALREADY_SET) err_strg='ELPA_ERROR_ENTRY_ALREADY_SET'
356  if (err_code_==ELPA_ERROR_ENTRY_INVALID_VALUE) err_strg='ELPA_ERROR_ENTRY_INVALID_VALUE'
357  if (err_code_==ELPA_ERROR_ENTRY_NO_STRING_REPRESENTATION) err_strg='ELPA_ERROR_NO_STRING_REPRESENTATION'
358 #endif
359 
360  write(msg,'(a)') 'ELPA library error!'
361  if (present(err_msg)) then
362    if (trim(err_msg)/="") write(msg,'(3a)') trim(msg),ch10,trim(err_msg)
363  end if
364  if (present(err_varname)) then
365    if (trim(err_varname)/="") write(msg,'(4a)') trim(msg),ch10,'Variable: ',trim(err_varname)
366  end if
367  if (trim(err_strg)/="") write(msg,'(4a)') trim(msg),ch10,'Error code: ',trim(err_strg)
368    ABI_ERROR(msg)
369 
370 end subroutine elpa_func_error_handler

m_elpa/elpa_func_get_communicators [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_get_communicators

FUNCTION

  Wrapper to elpa_get_communicators ELPA function

INPUTS

  mpi_comm_parent=Global communicator for the calculations (in)
  process_row=Row coordinate of the calling process in the process grid (in)
  process_col=Column coordinate of the calling process in the process grid (in)

SIDE EFFECTS

  elpa_hdl(type<elpa_hdl_t>)= ELPA handle

SOURCE

392 subroutine elpa_func_get_communicators(elpa_hdl,mpi_comm_parent,process_row,process_col)
393 
394 !Arguments ------------------------------------
395  integer,intent(in)  :: mpi_comm_parent,process_row,process_col
396  type(elpa_hdl_t),intent(inout) :: elpa_hdl
397 
398 !Local variables-------------------------------
399  integer  :: err
400  character(len=20) :: varname
401 
402 ! *********************************************************************
403 
404  err=0 ; varname=''
405 
406  if (.not.elpa_hdl%is_allocated) then
407    ABI_BUG('ELPA handle not allocated!')
408  end if
409 
410 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
411  if (err==ELPA_OK) then
412    varname='mpi_comm_parent'
413    call elpa_hdl%elpa%set(trim(varname),mpi_comm_parent,err)
414  end if
415  if (err==ELPA_OK) then
416    varname='process_row'
417    call elpa_hdl%elpa%set(trim(varname),process_row,err)
418  end if
419  if (err==ELPA_OK) then
420    varname='process_col'
421    call elpa_hdl%elpa%set(trim(varname),process_col,err)
422  end if
423  if (err==ELPA_OK) then
424    varname=''
425    err = elpa_hdl%elpa%setup()
426    call elpa_func_error_handler(err_code=err,err_msg='Error during ELPA setup')
427  endif
428 
429 #else
430  elpa_hdl%mpi_comm_parent=mpi_comm_parent
431  elpa_hdl%process_row=process_row
432  elpa_hdl%process_col=process_col
433 #if (defined HAVE_LINALG_ELPA_2016)
434  err=elpa_get_communicators(mpi_comm_parent,process_row,process_col,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
435 #elif (defined HAVE_LINALG_ELPA_2015_11) || (defined HAVE_LINALG_ELPA_2015_02)
436  err=get_elpa_row_col_comms(mpi_comm_parent,process_row,process_col,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
437 #elif (defined HAVE_LINALG_ELPA_2014) || (defined HAVE_LINALG_ELPA_2013)
438  call get_elpa_row_col_comms(mpi_comm_parent,process_row,process_col,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
439 #else
440 !ELPA-LEGACY-2017
441  err=elpa_get_communicators(mpi_comm_parent,process_row,process_col,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
442 #endif
443 
444 #endif
445 
446  call elpa_func_error_handler(err_code=err,err_msg='Error in elpa_get_communicators',err_varname=varname)
447 
448 end subroutine elpa_func_get_communicators

m_elpa/elpa_func_hermitian_multiply_complex [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_hermitian_multiply_complex

FUNCTION

  Wrapper to elpa_hermitian_multiply_complex ELPA function
  Performs C := A**H * B

INPUTS

 uplo_a='U' if A is upper triangular
        'L' if A is lower triangular
        anything else if A is a full matrix
        Please note: This pertains to the original A (as set in the calling program)
          whereas the transpose of A is used for calculations
          If uplo_a is 'U' or 'L', the other triangle is not used at all,
          i.e. it may contain arbitrary numbers
  uplo_c='U' if only the upper diagonal part of C is needed
         'L' if only the upper diagonal part of C is needed
         anything else if the full matrix C is needed
         Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
         written to a certain extent, i.e. one shouldn't rely on the content there!
  ncb=Number of columns of B and C
  aa(local_nrows,local_ncols)=Matrix A
  bb(ldb,local_ncols_c)=Matrix B
  local_nrows_b=Local rows of matrix B
  local_ncols_b=Local columns of matrix B
  local_nrows_c=Local rows of matrix C
  local_ncols_c=Local columns of matrix C

OUTPUT

  cc(local_nrows_c,local_ncols_c)=Matrix C

SIDE EFFECTS

  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

1297 subroutine elpa_func_hermitian_multiply_complex(elpa_hdl,uplo_a,uplo_c,ncb,aa,bb,local_nrows_b,local_ncols_b,&
1298 &                                               cc,local_nrows_c,local_ncols_c)
1299 
1300 !Arguments ------------------------------------
1301 !scalars
1302  integer,intent(in)  :: ncb,local_nrows_b,local_nrows_c,local_ncols_b,local_ncols_c
1303  character*1 :: uplo_a,uplo_c
1304  type(elpa_hdl_t),intent(inout) :: elpa_hdl
1305 !arrays
1306  complex(dpc),intent(in) :: aa(:,:),bb(:,:)
1307  complex(dpc),intent(out) :: cc(:,:)
1308 
1309 !Local variables-------------------------------
1310  integer :: err
1311  logical :: success
1312 
1313 ! *********************************************************************
1314 
1315  success=.true. ; err=0
1316 
1317  if (.not.elpa_hdl%is_allocated) then
1318    ABI_BUG('ELPA handle not allocated!')
1319  end if
1320  if (.not.elpa_hdl%matrix_is_set) then
1321    ABI_BUG('Matrix not set in ELPA handle!')
1322  end if
1323 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1324  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
1325 #else
1326  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
1327 #endif
1328  ABI_CHECK(size(bb)==local_nrows_b*local_ncols_b,'BUG: matrix B has wrong sizes!')
1329  ABI_CHECK(size(cc)==local_nrows_c*local_ncols_c,'BUG: matrix C has wrong sizes!')
1330 
1331 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1332  call elpa_hdl%elpa%hermitian_multiply(uplo_a,uplo_c,ncb,aa,bb,local_nrows_b,local_ncols_b,&
1333 &                                      cc,local_nrows_c,local_ncols_c,err)
1334  success=(err==ELPA_OK)
1335 #elif (defined HAVE_LINALG_ELPA_2016)
1336  success = elpa_mult_ah_b_complex(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,elpa_hdl%local_ncols,&
1337 &          bb,local_nrows_b,local_ncols_b,elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,&
1338 &          cc,local_nrows_c,local_ncols_c)
1339 #elif (defined HAVE_LINALG_ELPA_2015_11) || (defined HAVE_LINALG_ELPA_2015_02)
1340  call mult_ah_b_complex(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,bb,local_nrows_b,&
1341 &               elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,cc,local_nrows_c)
1342 #elif (defined HAVE_LINALG_ELPA_2014) || (defined HAVE_LINALG_ELPA_2013)
1343  call mult_ah_b_complex(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,bb,local_nrows_b,&
1344 &               elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,cc,local_nrows_c)
1345 #else
1346 !ELPA-LEGACY-2017
1347  success =  elpa_mult_ah_b_complex_double(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,&
1348 &           elpa_hdl%local_ncols,bb,local_nrows_b,local_ncols_b,elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,&
1349 &           elpa_hdl%elpa_comm_cols,cc,local_nrows_c,local_ncols_c)
1350 #endif
1351 
1352  if (.not.success) call elpa_func_error_handler(err_msg='Error in hermitian_multiply_complex!')
1353 
1354 end subroutine elpa_func_hermitian_multiply_complex

m_elpa/elpa_func_hermitian_multiply_real [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_hermitian_multiply_real

FUNCTION

  Wrapper to elpa_hermitian_multiply_real ELPA function
  Performs C := A**T * B

INPUTS

 uplo_a='U' if A is upper triangular
        'L' if A is lower triangular
        anything else if A is a full matrix
        Please note: This pertains to the original A (as set in the calling program)
          whereas the transpose of A is used for calculations
          If uplo_a is 'U' or 'L', the other triangle is not used at all,
          i.e. it may contain arbitrary numbers
  uplo_c='U' if only the upper diagonal part of C is needed
         'L' if only the upper diagonal part of C is needed
         anything else if the full matrix C is needed
         Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
         written to a certain extent, i.e. one shouldn't rely on the content there!
  ncb=Number of columns of B and C
  aa(local_nrows,local_ncols)=Matrix A
  bb(ldb,local_ncols_c)=Matrix B
  local_nrows_b=Local rows of matrix B
  local_ncols_b=Local columns of matrix B
  local_nrows_c=Local rows of matrix C
  local_ncols_c=Local columns of matrix C

OUTPUT

  cc(local_nrows_c,local_ncols_c)=Matrix C

SIDE EFFECTS

  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

1198 subroutine elpa_func_hermitian_multiply_real(elpa_hdl,uplo_a,uplo_c,ncb,aa,bb,local_nrows_b,local_ncols_b,&
1199 &                                            cc,local_nrows_c,local_ncols_c)
1200 
1201 !Arguments ------------------------------------
1202 !scalars
1203  integer,intent(in)  :: ncb,local_nrows_b,local_nrows_c,local_ncols_b,local_ncols_c
1204  character*1 :: uplo_a,uplo_c
1205  type(elpa_hdl_t),intent(inout) :: elpa_hdl
1206 !arrays
1207  real(dp),intent(in) :: aa(:,:),bb(:,:)
1208  real(dp),intent(out) :: cc(:,:)
1209 
1210 !Local variables-------------------------------
1211  integer :: err
1212  logical :: success
1213 
1214 ! *********************************************************************
1215 
1216  success=.true. ; err=0
1217 
1218  if (.not.elpa_hdl%is_allocated) then
1219    ABI_BUG('ELPA handle not allocated!')
1220  end if
1221  if (.not.elpa_hdl%matrix_is_set) then
1222    ABI_BUG('Matrix not set in ELPA handle!')
1223  end if
1224 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1225  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
1226 #else
1227  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
1228 #endif
1229  ABI_CHECK(size(bb)==local_nrows_b*local_ncols_b,'BUG: matrix B has wrong sizes!')
1230  ABI_CHECK(size(cc)==local_nrows_c*local_ncols_c,'BUG: matrix C has wrong sizes!')
1231 
1232 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1233  call elpa_hdl%elpa%hermitian_multiply(uplo_a,uplo_c,ncb,aa,bb,local_nrows_b,local_ncols_b,&
1234 &                                      cc,local_nrows_c,local_ncols_c,err)
1235  success=(err==ELPA_OK)
1236 #elif (defined HAVE_LINALG_ELPA_2016)
1237  success = elpa_mult_at_b_real(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,elpa_hdl%local_ncols,&
1238 &          bb,local_nrows_b,local_ncols_b,elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,&
1239 &          cc,local_nrows_c,local_ncols_c)
1240 #elif (defined HAVE_LINALG_ELPA_2015_11) || (defined HAVE_LINALG_ELPA_2015_02)
1241  call mult_at_b_real(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,bb,local_nrows_b,&
1242 &               elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,cc,local_nrows_c)
1243 #elif (defined HAVE_LINALG_ELPA_2014) || (defined HAVE_LINALG_ELPA_2013)
1244  call mult_at_b_real(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,bb,local_nrows_b,&
1245 &               elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,cc,local_nrows_c)
1246 #else
1247 !ELPA-LEGACY-2017
1248  success =  elpa_mult_at_b_real_double(uplo_a,uplo_c,elpa_hdl%na,ncb,aa,elpa_hdl%local_nrows,&
1249 &           elpa_hdl%local_ncols,bb,local_nrows_b,local_ncols_b,elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,&
1250 &           elpa_hdl%elpa_comm_cols,cc,local_nrows_c,local_ncols_c)
1251 #endif
1252 
1253  if (.not.success) call elpa_func_error_handler(err_msg='Error in hermitian_multiply_real!')
1254 
1255 end subroutine elpa_func_hermitian_multiply_real

m_elpa/elpa_func_init [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_init

FUNCTION

  Wrapper to elpa_init ELPA function

INPUTS

OUTPUT

SOURCE

126 subroutine elpa_func_init()
127 
128 !Local variables-------------------------------
129 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
130  integer,parameter :: elpa_min_version=20170403
131 #endif
132  logical :: success
133 !
134 ! *********************************************************************
135 
136  success=.true.
137 
138 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
139 
140 #if (defined HAVE_LINALG_ELPA_2016) || (defined HAVE_LINALG_ELPA_2015_11) || (defined HAVE_LINALG_ELPA_2015_02)
141 !This case is not supposed to happen
142  success=.false.
143  ABI_BUG('Wrong ELPA cpp directives!')
144 #elif (defined HAVE_LINALG_ELPA_2014) || (defined HAVE_LINALG_ELPA_2013)
145 !This case is not supposed to happen
146  success=.false.
147  ABI_BUG('Wrong ELPA cpp directives!')
148 #else
149  success=(elpa_init(elpa_min_version)==ELPA_OK)
150 #endif
151 
152 #else
153 !No init function
154 #endif
155 
156  if (.not.success) then
157    ABI_ERROR('Problem with ELPA (elpa_init)!')
158  end if
159 
160 end subroutine elpa_func_init

m_elpa/elpa_func_invert_triangular_complex [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_invert_triangular_complex

FUNCTION

  Wrapper to elpa_invert_triangular_complex ELPA function

INPUTS

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix which should be factorized.
                     Distribution is like in Scalapack.
                     Only upper triangle is needs to be set.
                     The lower triangle is not referenced.
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

1101 subroutine elpa_func_invert_triangular_complex(elpa_hdl,aa)
1102 
1103 !Arguments ------------------------------------
1104 !scalars
1105  type(elpa_hdl_t),intent(inout) :: elpa_hdl
1106 !arrays
1107  complex(dpc),intent(inout) :: aa(:,:)
1108 
1109 !Local variables-------------------------------
1110  integer :: err
1111  logical :: success
1112 
1113 ! *********************************************************************
1114 
1115  success=.true. ; err=0
1116 
1117  if (.not.elpa_hdl%is_allocated) then
1118    ABI_BUG('ELPA handle not allocated!')
1119  end if
1120  if (.not.elpa_hdl%matrix_is_set) then
1121    ABI_BUG('Matrix not set in ELPA handle!')
1122  end if
1123 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1124  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
1125 #else
1126  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
1127 #endif
1128 
1129 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1130  call elpa_hdl%elpa%invert_triangular(aa,err)
1131  success=(err==ELPA_OK)
1132 #elif (defined HAVE_LINALG_ELPA_2016)
1133  success = elpa_invert_trm_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
1134 &                                  elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
1135 #elif (defined HAVE_LINALG_ELPA_2015_11)
1136  call invert_trm_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
1137                          elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
1138 #elif (defined HAVE_LINALG_ELPA_2015_02)
1139  call invert_trm_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
1140                          elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
1141 #elif (defined HAVE_LINALG_ELPA_2014)
1142  call invert_trm_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
1143 &                        elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,success)
1144  success=.true. ! Sometimes get unexpected success=false
1145 #elif (defined HAVE_LINALG_ELPA_2013)
1146  call invert_trm_complex(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
1147 &                        elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
1148 #else
1149 !ELPA-LEGACY-2017
1150  success = elpa_invert_trm_complex_double(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
1151 &                                         elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
1152 #endif
1153 
1154  if (.not.success) call elpa_func_error_handler(err_msg='Error in invert_trianguler_complex!')
1155 
1156 end subroutine elpa_func_invert_triangular_complex

m_elpa/elpa_func_invert_triangular_real [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_invert_triangular_real

FUNCTION

  Wrapper to elpa_invert_triangular_real ELPA function

INPUTS

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix which should be factorized.
                     Distribution is like in Scalapack.
                     Only upper triangle is needs to be set.
                     The lower triangle is not referenced.
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

1023 subroutine elpa_func_invert_triangular_real(elpa_hdl,aa)
1024 
1025 !Arguments ------------------------------------
1026 !scalars
1027  type(elpa_hdl_t),intent(inout) :: elpa_hdl
1028 !arrays
1029  real(dp),intent(inout) :: aa(:,:)
1030 
1031 !Local variables-------------------------------
1032  integer :: err
1033  logical :: success
1034 
1035 ! *********************************************************************
1036 
1037  success=.true. ; err=0
1038 
1039  if (.not.elpa_hdl%is_allocated) then
1040    ABI_BUG('ELPA handle not allocated!')
1041  end if
1042  if (.not.elpa_hdl%matrix_is_set) then
1043    ABI_BUG('Matrix not set in ELPA handle!')
1044  end if
1045 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1046  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
1047 #else
1048  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
1049 #endif
1050 
1051 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
1052  call elpa_hdl%elpa%invert_triangular(aa,err)
1053  success=(err==ELPA_OK)
1054 #elif (defined HAVE_LINALG_ELPA_2016)
1055  success = elpa_invert_trm_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
1056 &                               elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
1057 #elif (defined HAVE_LINALG_ELPA_2015_11)
1058  call invert_trm_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
1059                          elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
1060 #elif (defined HAVE_LINALG_ELPA_2015_02)
1061  call invert_trm_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
1062                          elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.,success)
1063 #elif (defined HAVE_LINALG_ELPA_2014)
1064  call invert_trm_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
1065 &                     elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,success)
1066  success=.true. ! Sometimes get unexpected success=false
1067 #elif (defined HAVE_LINALG_ELPA_2013)
1068  call invert_trm_real(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,&
1069 &                     elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
1070 #else
1071 !ELPA-LEGACY-2017
1072  success = elpa_invert_trm_real_double(elpa_hdl%na,aa,elpa_hdl%local_nrows,elpa_hdl%nblk,elpa_hdl%local_ncols,&
1073 &                                      elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,.false.)
1074 #endif
1075 
1076  if (.not.success) call elpa_func_error_handler(err_msg='Error in invert_trianguler_real!')
1077 
1078 end subroutine elpa_func_invert_triangular_real

m_elpa/elpa_func_set_matrix [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_set_matrix

FUNCTION

  Set parameters decribing a matrix and it's MPI distribution
  in a ELPA handle

INPUTS

  na=Order of matrix A
  nblk=Blocksize of cyclic distribution, must be the same in both directions!
  nev=Number of eigenvalues needed.
  local_nrows=Leading dimension of A
  local_ncols=Local columns of matrixes A and Q (eigenvectors)
  nev=Number of eigenvalues needed.

SIDE EFFECTS

  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

474 subroutine elpa_func_set_matrix(elpa_hdl,na,nblk,nev,local_nrows,local_ncols)
475 
476 !Arguments ------------------------------------
477 !scalars
478  integer,intent(in) :: na,nblk,nev,local_nrows,local_ncols
479  type(elpa_hdl_t),intent(inout) :: elpa_hdl
480 !arrays
481 
482 !Local variables-------------------------------
483  integer :: err
484  character(len=15) :: varname
485 
486 ! *********************************************************************
487 
488  err=0 ; varname=''
489 
490  if (.not.elpa_hdl%is_allocated) then
491    ABI_BUG('ELPA handle not allocated!')
492  end if
493 
494 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
495  if (err==ELPA_OK) then
496    varname="na"
497    call elpa_hdl%elpa%set(trim(varname),na,err)
498  end if
499  if (err==ELPA_OK) then
500    varname="nblk"
501    call elpa_hdl%elpa%set(trim(varname),nblk,err)
502  end if
503  if (err==ELPA_OK) then
504    varname="nev"
505    call elpa_hdl%elpa%set(trim(varname),nev,err)
506  end if
507  if (err==ELPA_OK) then
508    varname="local_nrows"
509    call elpa_hdl%elpa%set(trim(varname),local_nrows,err)
510  end if
511  if (err==ELPA_OK) then
512    varname="local_ncols"
513    call elpa_hdl%elpa%set(trim(varname),local_ncols,err)
514  end if
515 #else
516  elpa_hdl%na=na
517  elpa_hdl%nblk=nblk
518  elpa_hdl%nev=nev
519  elpa_hdl%local_nrows=local_nrows
520  elpa_hdl%local_ncols=local_ncols
521 #endif
522 
523  call elpa_func_error_handler(err_code=err,err_msg='Error during matrix initialization',err_varname=varname)
524 
525  elpa_hdl%matrix_is_set=.true.
526 
527 end subroutine elpa_func_set_matrix

m_elpa/elpa_func_solve_evp_1stage_complex [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_solve_evp_1stage_complex

FUNCTION

  Wrapper to elpa_solve_evp_complex_1stage ELPA function

INPUTS

  nev=Number of eigenvalues needed.

OUTPUT

  ev(na)=Eigenvalues of a, every processor gets the complete set
  qq(local_nrows,local_ncols)=Eigenvectors of aa
                     Distribution is like in Scalapack.
                     Must be always dimensioned to the full size (corresponding to (na,na))
                      even if only a part of the eigenvalues is needed.

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix for which eigenvalues are to be computed.
                    Distribution is like in Scalapack.
                    The full matrix must be set (not only one half like in scalapack).
                    Destroyed on exit (upper and lower half).
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

786 subroutine elpa_func_solve_evp_1stage_complex(elpa_hdl,aa,qq,ev,nev)
787 
788 !Arguments ------------------------------------
789 !scalars
790  integer,intent(in)  :: nev
791  type(elpa_hdl_t),intent(inout) :: elpa_hdl
792 !arrays
793  complex(dpc),intent(inout) :: aa(:,:)
794  real(dp),intent(out) :: ev(:)
795  complex(dpc),intent(out) :: qq(:,:)
796 
797 !Local variables-------------------------------
798  integer :: err
799  logical  :: success
800 
801 ! *********************************************************************
802 
803  success=.true. ; err=0
804 
805  if (.not.elpa_hdl%is_allocated) then
806    ABI_BUG('ELPA handle not allocated!')
807  end if
808  if (.not.elpa_hdl%matrix_is_set) then
809    ABI_BUG('Matrix not set in ELPA handle!')
810  end if
811 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
812  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
813  ABI_CHECK(size(qq)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix Q has wrong sizes!')
814  ABI_CHECK(size(ev)==elpa_hdl%elpa%na,'BUG: matrix EV has wrong sizes!')
815 #else
816  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
817  ABI_CHECK(size(qq)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix Q has wrong sizes!')
818  ABI_CHECK(size(ev)==elpa_hdl%na,'BUG: matrix EV has wrong sizes!')
819 #endif
820 
821 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
822  if (err==ELPA_OK) call elpa_hdl%elpa%set("solver",ELPA_SOLVER_1STAGE,err)
823  if (err==ELPA_OK) call elpa_hdl%elpa%eigenvectors(aa,ev,qq,err)
824  success=(err==ELPA_OK)
825 #elif  (defined HAVE_LINALG_ELPA_2016)
826  success=elpa_solve_evp_complex_1stage(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
827 &                       elpa_hdl%nblk,elpa_hdl%local_ncols,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
828 #elif (defined HAVE_LINALG_ELPA_2015_11)
829  success=solve_evp_complex(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
830 &                       elpa_hdl%nblk,elpa_hdl%local_ncols,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
831 #elif (defined HAVE_LINALG_ELPA_2015_02) || (defined HAVE_LINALG_ELPA_2014)
832  success=solve_evp_complex(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
833 &                       elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
834 #elif (defined HAVE_LINALG_ELPA_2013)
835  call solve_evp_complex(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
836 &                    elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
837 #else
838 !ELPA-LEGACY-2017
839  success=elpa_solve_evp_complex_1stage_double(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
840 &        elpa_hdl%nblk,elpa_hdl%local_ncols,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,elpa_hdl%mpi_comm_parent,.false.)
841 #endif
842 
843  if (.not.success) call elpa_func_error_handler(err_msg='Error in solve_evp_1stage_complex!')
844 
845 end subroutine elpa_func_solve_evp_1stage_complex

m_elpa/elpa_func_solve_evp_1stage_real [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_solve_evp_1stage_real

FUNCTION

  Wrapper to elpa_solve_evp_real_1stage ELPA function

INPUTS

  nev=Number of eigenvalues needed.

OUTPUT

  ev(na)=Eigenvalues of a, every processor gets the complete set
  qq(local_nrows,local_ncols)=Eigenvectors of aa
                     Distribution is like in Scalapack.
                     Must be always dimensioned to the full size (corresponding to (na,na))
                     even if only a part of the eigenvalues is needed.

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix for which eigenvalues are to be computed.
                    Distribution is like in Scalapack.
                    The full matrix must be set (not only one half like in scalapack).
                    Destroyed on exit (upper and lower half).
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

697 subroutine elpa_func_solve_evp_1stage_real(elpa_hdl,aa,qq,ev,nev)
698 
699 !Arguments ------------------------------------
700 !scalars
701  integer,intent(in)  :: nev
702  type(elpa_hdl_t),intent(inout) :: elpa_hdl
703 !arrays
704  real(dp),intent(inout) :: aa(:,:)
705  real(dp),intent(out) :: ev(:),qq(:,:)
706 
707 !Local variables-------------------------------
708  integer :: err
709  logical  :: success
710 
711 ! *********************************************************************
712 
713  success=.true. ; err=0
714 
715  if (.not.elpa_hdl%is_allocated) then
716    ABI_BUG('ELPA handle not allocated!')
717  end if
718  if (.not.elpa_hdl%matrix_is_set) then
719    ABI_BUG('Matrix not set in ELPA handle!')
720  end if
721 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
722  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
723  ABI_CHECK(size(qq)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix Q has wrong sizes!')
724  ABI_CHECK(size(ev)==elpa_hdl%elpa%na,'BUG: matrix EV has wrong sizes!')
725 #else
726  ABI_CHECK(size(aa)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix A has wrong sizes!')
727  ABI_CHECK(size(qq)==elpa_hdl%local_nrows*elpa_hdl%local_ncols,'BUG: matrix Q has wrong sizes!')
728  ABI_CHECK(size(ev)==elpa_hdl%na,'BUG: matrix EV has wrong sizes!')
729 #endif
730 
731 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
732  if (err==ELPA_OK) call elpa_hdl%elpa%set("solver",ELPA_SOLVER_1STAGE,err)
733  if (err==ELPA_OK) call elpa_hdl%elpa%eigenvectors(aa,ev,qq,err)
734  success=(err==ELPA_OK)
735 #elif  (defined HAVE_LINALG_ELPA_2016)
736  success=elpa_solve_evp_real_1stage(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
737 &                       elpa_hdl%nblk,elpa_hdl%local_ncols,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
738 #elif (defined HAVE_LINALG_ELPA_2015_11)
739  success=solve_evp_real(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
740 &                       elpa_hdl%nblk,elpa_hdl%local_ncols,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
741 #elif (defined HAVE_LINALG_ELPA_2015_02) || (defined HAVE_LINALG_ELPA_2014)
742  success=solve_evp_real(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
743 &                       elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
744 #elif (defined HAVE_LINALG_ELPA_2013)
745  call solve_evp_real(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
746 &                    elpa_hdl%nblk,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols)
747 #else
748 !ELPA-LEGACY-2017
749  success=elpa_solve_evp_real_1stage_double(elpa_hdl%na,nev,aa,elpa_hdl%local_nrows,ev,qq,elpa_hdl%local_nrows,&
750 &        elpa_hdl%nblk,elpa_hdl%local_ncols,elpa_hdl%elpa_comm_rows,elpa_hdl%elpa_comm_cols,elpa_hdl%mpi_comm_parent,.false.)
751 #endif
752 
753  if (.not.success) call elpa_func_error_handler(err_msg='Error in solve_evp_1stage_real!')
754 
755 end subroutine elpa_func_solve_evp_1stage_real

m_elpa/elpa_func_solve_evp_2stage_complex [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_solve_evp_2stage_complex

FUNCTION

  Wrapper to elpa_solve_evp_complex_2stage ELPA function

INPUTS

  nev=Number of eigenvalues needed.

OUTPUT

  ev(na)=Eigenvalues of a, every processor gets the complete set
  qq(local_nrows,local_ncols)=Eigenvectors of aa
                     Distribution is like in Scalapack.
                     Must be always dimensioned to the full size (corresponding to (na,na))
                      even if only a part of the eigenvalues is needed.

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix for which eigenvalues are to be computed.
                    Distribution is like in Scalapack.
                    The full matrix must be set (not only one half like in scalapack).
                    Destroyed on exit (upper and lower half).
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

627 subroutine elpa_func_solve_gevp_2stage_complex(elpa_hdl,aa,bb,qq,ev,nev)
628 
629 !Arguments ------------------------------------
630 !scalars
631  integer,intent(in)  :: nev
632  type(elpa_hdl_t),intent(inout) :: elpa_hdl
633 !arrays
634  complex(dpc),intent(inout) :: aa(:,:),bb(:,:)
635  real(dp),intent(out) :: ev(:)
636  complex(dpc),intent(out) :: qq(:,:)
637 
638 !Local variables-------------------------------
639  integer :: err
640  logical  :: success
641 
642 ! *********************************************************************
643 
644  success=.true. ; err=0
645 
646  if (.not.elpa_hdl%is_allocated) then
647    ABI_BUG('ELPA handle not allocated!')
648  end if
649  if (.not.elpa_hdl%matrix_is_set) then
650    ABI_BUG('Matrix not set in ELPA handle!')
651  end if
652 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
653  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
654  ABI_CHECK(size(qq)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix Q has wrong sizes!')
655  ABI_CHECK(size(ev)==elpa_hdl%elpa%na,'BUG: matrix EV has wrong sizes!')
656 #endif
657 
658 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
659  if (err==ELPA_OK) call elpa_hdl%elpa%set("solver",ELPA_SOLVER_2STAGE,err)
660  if (err==ELPA_OK) call elpa_hdl%elpa%generalized_eigenvectors(aa,bb,ev,qq,.false.,err)
661  success=(err==ELPA_OK)
662 #endif
663 
664  if (.not.success) call elpa_func_error_handler(err_code=err,err_msg='Error in solve_evp_2stage_complex!')
665 
666 end subroutine elpa_func_solve_gevp_2stage_complex

m_elpa/elpa_func_solve_evp_2stage_real [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_solve_evp_2stage_real

FUNCTION

  Wrapper to elpa_solve_evp_real_2stage ELPA function

INPUTS

  nev=Number of eigenvalues needed.

OUTPUT

  ev(na)=Eigenvalues of a, every processor gets the complete set
  qq(local_nrows,local_ncols)=Eigenvectors of aa
                     Distribution is like in Scalapack.
                     Must be always dimensioned to the full size (corresponding to (na,na))
                     even if only a part of the eigenvalues is needed.

SIDE EFFECTS

  aa(local_nrows,local_ncols)=Distributed matrix for which eigenvalues are to be computed.
                    Distribution is like in Scalapack.
                    The full matrix must be set (not only one half like in scalapack).
                    Destroyed on exit (upper and lower half).
  elpa_hdl(type<elpa_hdl_t>)=handler for ELPA object

SOURCE

558 subroutine elpa_func_solve_gevp_2stage_real(elpa_hdl,aa,bb,qq,ev,nev)
559 
560 !Arguments ------------------------------------
561 !scalars
562  integer,intent(in)  :: nev
563  type(elpa_hdl_t),intent(inout) :: elpa_hdl
564 !arrays
565  real(dp),intent(inout) :: aa(:,:),bb(:,:)
566  real(dp),intent(out) :: ev(:),qq(:,:)
567 
568 !Local variables-------------------------------
569  integer :: err
570  logical  :: success
571 
572 ! *********************************************************************
573 
574  success=.true. ; err=0
575 
576  if (.not.elpa_hdl%is_allocated) then
577    ABI_BUG('ELPA handle not allocated!')
578  end if
579  if (.not.elpa_hdl%matrix_is_set) then
580    ABI_BUG('Matrix not set in ELPA handle!')
581  end if
582 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
583  ABI_CHECK(size(aa)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix A has wrong sizes!')
584  ABI_CHECK(size(qq)==elpa_hdl%elpa%local_nrows*elpa_hdl%elpa%local_ncols,'BUG: matrix Q has wrong sizes!')
585  ABI_CHECK(size(ev)==elpa_hdl%elpa%na,'BUG: matrix EV has wrong sizes!')
586 #endif
587 
588 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
589  if (err==ELPA_OK) call elpa_hdl%elpa%set("solver",ELPA_SOLVER_2STAGE,err)
590  if (err==ELPA_OK) call elpa_hdl%elpa%generalized_eigenvectors(aa,bb,ev,qq,.false.,err)
591  success=(err==ELPA_OK)
592 #endif
593 
594  if (.not.success) call elpa_func_error_handler(err_msg='Error in solve_evp_2stage_real!')
595 
596 end subroutine elpa_func_solve_gevp_2stage_real

m_elpa/elpa_func_uninit [ Functions ]

[ Top ] [ m_elpa ] [ Functions ]

NAME

  elpa_func_uninit

FUNCTION

  Wrapper to elpa_uninit ELPA function

INPUTS

OUTPUT

SOURCE

178 subroutine elpa_func_uninit()
179 
180 ! *********************************************************************
181 
182 #ifdef HAVE_LINALG_ELPA_FORTRAN2008
183  call elpa_uninit()
184 #else
185 !No uninit function
186 #endif
187 
188 end subroutine elpa_func_uninit