Line data Source code
1 : !-------------------------------------------------------------------------------
2 : ! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
3 : ! This file is part of FLEUR and available as free software under the conditions
4 : ! of the MIT license as expressed in the LICENSE file in more detail.
5 : !--------------------------------------------------------------------------------
6 :
7 : MODULE m_elpa_onenode
8 : CONTAINS
9 0 : SUBROUTINE elpa_diag_onenode(hmat, smat, ne, eig, ev)
10 : !
11 : !----------------------------------------------------
12 : ! Sigensystem solver
13 : ! Uses the ELPA (1 node version)
14 : !
15 : !
16 : ! hmat ..... Hamiltonian matrix
17 : ! smat ..... overlap matrix
18 : ! ne ....... number of ev's searched (and found) on this node
19 : ! On input, overall number of ev's searched,
20 : ! On output, local number of ev's found
21 : ! eig ...... eigenvalues, output
22 : ! ev ....... eigenvectors, output
23 : !
24 : ! U.Alekseeva Nov. 2018
25 : !----------------------------------------------------
26 : USE m_juDFT
27 : USE m_types_mat
28 : #ifdef CPP_ELPA_ONENODE
29 : USE elpa
30 :
31 : #endif
32 : IMPLICIT NONE
33 :
34 : CLASS(t_mat), INTENT(INOUT) :: hmat, smat
35 : CLASS(t_mat), ALLOCATABLE, INTENT(OUT)::ev
36 : REAL, INTENT(OUT) :: eig(:)
37 : INTEGER, INTENT(INOUT) :: ne
38 :
39 : #ifdef CPP_ELPA_ONENODE
40 : !... Local variables
41 : !
42 : INTEGER :: err
43 : REAL, ALLOCATABLE :: eig2(:)
44 : TYPE(t_mat) :: ev_dist
45 : INTEGER :: kernel
46 : CLASS(elpa_t), pointer :: elpa_obj
47 :
48 : call timestart("ELPA 2018 one-node")
49 :
50 :
51 : err = elpa_init(20180525)
52 : elpa_obj => elpa_allocate()
53 :
54 : ALLOCATE (eig2(hmat%matsize1), stat=err) ! The eigenvalue array
55 : IF (err .NE. 0) CALL juDFT_error('Failed to allocated "eig2"', calledby='elpa')
56 :
57 : CALL ev_dist%init(hmat)! Eigenvectors
58 : IF (err .NE. 0) CALL juDFT_error('Failed to allocated "ev_dist"', calledby='elpa')
59 :
60 : CALL elpa_obj%set("na", hmat%matsize1, err)
61 : CALL elpa_obj%set("nev", ne, err)
62 : CALL elpa_obj%set("local_nrows", hmat%matsize1, err)
63 : CALL elpa_obj%set("local_ncols", hmat%matsize2, err)
64 : CALL elpa_obj%set("nblk", hmat%matsize1, err)
65 : CALL elpa_obj%set("blacs_context", -1, err)
66 : #ifdef CPP_GPU
67 : CALL elpa_obj%set("gpu", 1, err)
68 : #endif
69 : err = elpa_obj%setup()
70 : call elpa_obj%get("solver", kernel)
71 : print *, "elpa uses "//elpa_int_value_to_string("solver", kernel)//" solver"
72 :
73 : CALL hmat%add_transpose(hmat)
74 : CALL smat%add_transpose(smat)
75 :
76 : IF (hmat%l_real) THEN
77 : CALL elpa_obj%generalized_eigenvectors(hmat%data_r, smat%data_r, eig2, ev_dist%data_r, .FALSE., err)
78 : ELSE
79 : CALL elpa_obj%generalized_eigenvectors(hmat%data_c, smat%data_c, eig2, ev_dist%data_c, .FALSE., err)
80 : ENDIF
81 :
82 :
83 : CALL elpa_deallocate(elpa_obj)
84 : CALL elpa_uninit()
85 : ! END of ELPA stuff
86 :
87 : eig(1:ne) = eig2(1:ne)
88 : DEALLOCATE (eig2)
89 :
90 : ALLOCATE (t_mat::ev)
91 : CALL ev%alloc(hmat%l_real, hmat%matsize1, ne)
92 : CALL ev%copy(ev_dist, 1, 1)
93 :
94 :
95 : #endif
96 0 : call timestop("ELPA 2018 one-node")
97 :
98 0 : END SUBROUTINE elpa_diag_onenode
99 : END MODULE m_elpa_onenode
|