Line data Source code
1 : ! Copyright (c) 2018 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
2 : ! This file is part of FLEUR and available as free software under the conditions
3 : ! of the MIT license as expressed in the LICENSE file in more detail.
4 : !
5 : ! @authors: Miriam Hinzen, Gregor Michalicek
6 : ! Added MPI implementation, DW 2018
7 : !--------------------------------------------------------------------------------
8 : MODULE m_dummy_diag
9 : USE m_judft
10 : USE m_constants
11 : IMPLICIT NONE
12 : PRIVATE
13 :
14 : PUBLIC dummy_diag
15 :
16 : CONTAINS
17 :
18 0 : SUBROUTINE dummy_diag(hmat,smat,ne,eig,zmat)
19 : !Dummy diver: does not solve actual eigenvalue problem but simply returns a set of orthogonal vectors.
20 : !Could be useful for performance testing workloads in which we do not want to look at the diagonalization.
21 : ! A Cholesky decomp is still done to be able to do a back transform so that the resulting vector are orthonormal
22 : ! with respect to overlapp matrix.
23 :
24 : USE m_types_mat
25 : USE m_judft
26 :
27 : IMPLICIT NONE
28 :
29 : TYPE(t_mat), INTENT(INOUT) :: hmat,smat
30 : INTEGER, INTENT(INOUT) :: ne
31 : CLASS(t_mat), ALLOCATABLE, INTENT(OUT) :: zmat
32 : REAL, INTENT(OUT) :: eig(:)
33 :
34 : INTEGER :: nev,lwork,liwork,n
35 : INTEGER :: info
36 :
37 :
38 :
39 0 : ALLOCATE(t_mat::zmat)
40 0 : CALL zmat%alloc(hmat%l_real,hmat%matsize1,ne)
41 :
42 :
43 0 : IF (hmat%l_real) THEN
44 : ! --> start with Cholesky factorization of b ( so that b = l * l^t)
45 : ! --> b is overwritten by l
46 0 : CALL dpotrf('U',smat%matsize1,smat%data_r,SIZE(smat%data_r,1),info)
47 0 : IF (info.NE.0) THEN
48 0 : WRITE (*,*) 'Error in dpotrf: info =',info
49 0 : CALL juDFT_error("Diagonalization failed",calledby="lapack_singlePrec_diag")
50 : ENDIF
51 :
52 :
53 : ! --> solve a' * z' = eig * z' for eigenvalues eig between lb und ub
54 0 : zmat%data_r=0.0
55 0 : DO n=1,ne
56 0 : eig(ne)=-0.1+ne*1E-5
57 0 : zmat%data_r(ne,ne)=1.0
58 : enddo
59 : ! --> recover the generalized eigenvectors z by solving z' = l^t * z
60 0 : CALL dtrtrs('U','N','N',hmat%matsize1,nev,smat%data_r,smat%matsize1,zMat%data_r,zmat%matsize1,info)
61 0 : IF (info.NE.0) THEN
62 0 : WRITE (oUnit,*) 'Error in dtrtrs: info =',info
63 0 : CALL juDFT_error("Diagonalization failed",calledby="lapack_singlePrec_diag")
64 : ENDIF
65 :
66 :
67 : ELSE
68 :
69 : ! --> start with Cholesky factorization of b ( so that b = l * l^t)
70 : ! --> b is overwritten by l
71 0 : CALL zpotrf('U',smat%matsize1,smat%data_c,SIZE(smat%data_c,1),info)
72 0 : IF (info.NE.0) THEN
73 0 : WRITE (*,*) 'Error in zpotrf: info =',info
74 0 : CALL juDFT_error("Diagonalization failed",calledby="chase_diag")
75 : ENDIF
76 :
77 :
78 : ! --> solve a' * z' = eig * z' for eigenvalues eig between lb und ub
79 0 : zmat%data_c=0.0
80 0 : DO n=1,ne
81 0 : eig(ne)=-0.1+ne*1E-5
82 0 : zmat%data_c(ne,ne)=1.0
83 : enddo
84 :
85 : ! --> recover the generalized eigenvectors z by solving z' = l^t * z
86 0 : CALL ztrtrs('U','N','N',hmat%matsize1,nev,smat%data_c,smat%matsize1,zMat%data_c,zmat%matsize1,info)
87 0 : IF (info.NE.0) THEN
88 0 : WRITE (oUnit,*) 'Error in ztrtrs: info =',info
89 0 : CALL juDFT_error("Diagonalization failed",calledby="chase_diag")
90 : ENDIF
91 :
92 :
93 : ENDIF
94 0 : END SUBROUTINE dummy_diag
95 :
96 : END MODULE m_dummy_diag
|