(git:d18deda)
Loading...
Searching...
No Matches
cp_cfm_dlaf_api.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
9
11 USE cp_cfm_types, ONLY: cp_cfm_type
12#if defined(__DLAF)
13 USE dlaf_fortran, ONLY: dlaf_pzheevd, &
14 dlaf_pzhegvd, &
15 dlaf_pzpotrf
16#endif
17 USE kinds, ONLY: dp
18#include "../base/base_uses.f90"
19
20 IMPLICIT NONE
21
22 PRIVATE
23
24 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp_cfm_dlaf_api'
25
26 PUBLIC :: cp_cfm_pzpotrf_dlaf
28
29CONTAINS
30
31!***************************************************************************************************
32!> \brief Cholesky factorization using DLA-Future
33!> \param uplo ...
34!> \param n Matrix size
35!> \param a Local matrix
36!> \param ia Row index of first row (has to be 1)
37!> \param ja Col index of first column (has to be 1)
38!> \param desca ScaLAPACK matrix descriptor
39!> \param info 0 if factorization completed normally
40!> \author Rocco Meli
41! **************************************************************************************************
42 SUBROUTINE cp_cfm_pzpotrf_dlaf(uplo, n, a, ia, ja, desca, info)
43 CHARACTER, INTENT(IN) :: uplo
44 INTEGER, INTENT(IN) :: n
45 COMPLEX(KIND=dp), DIMENSION(:, :), TARGET :: a
46 INTEGER, INTENT(IN) :: ia, ja
47 INTEGER, DIMENSION(9) :: desca
48 INTEGER, TARGET :: info
49
50 CHARACTER(len=*), PARAMETER :: routinen = 'cp_cfm_pzpotrf_dlaf'
51
52 INTEGER :: handle
53
54 CALL timeset(routinen, handle)
55#if defined(__DLAF)
56 CALL dlaf_pzpotrf(uplo, n, a, ia, ja, desca, info)
57#else
58 mark_used(uplo)
59 mark_used(n)
60 mark_used(a)
61 mark_used(ia)
62 mark_used(ja)
63 mark_used(desca)
64 mark_used(info)
65 cpabort("CP2K compiled without the DLA-Future library.")
66#endif
67 CALL timestop(handle)
68 END SUBROUTINE cp_cfm_pzpotrf_dlaf
69
70 ! **************************************************************************************************
71!> \brief DLA-Future eigensolver for complex Hermitian matrices
72!> \param matrix ...
73!> \param eigenvectors ...
74!> \param eigenvalues ...
75!> \author Rocco Meli
76! **************************************************************************************************
77 SUBROUTINE cp_cfm_diag_dlaf(matrix, eigenvectors, eigenvalues)
78
79 TYPE(cp_cfm_type), INTENT(IN) :: matrix, eigenvectors
80 REAL(kind=dp), DIMENSION(:), INTENT(OUT) :: eigenvalues
81
82 CHARACTER(LEN=*), PARAMETER :: routinen = 'cp_cfm_diag_dlaf'
83
84 INTEGER :: handle, n, nmo
85 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), TARGET :: eig
86
87 CALL timeset(routinen, handle)
88
89 n = matrix%matrix_struct%nrow_global
90 ALLOCATE (eig(n))
91
92 CALL cp_cfm_diag_dlaf_base(matrix, eigenvectors, eig)
93
94 nmo = SIZE(eigenvalues, 1)
95 IF (nmo > n) THEN
96 eigenvalues(1:n) = eig(1:n)
97 ELSE
98 eigenvalues(1:nmo) = eig(1:nmo)
99 END IF
100
101 DEALLOCATE (eig)
102
103 CALL timestop(handle)
104
105 END SUBROUTINE cp_cfm_diag_dlaf
106
107! **************************************************************************************************
108!> \brief DLA-Future generalized eigensolver for complex Hermitian matrices
109!> \param amatrix ...
110!> \param bmatrix ...
111!> \param eigenvectors ...
112!> \param eigenvalues ...
113!> \author Rocco Meli
114! **************************************************************************************************
115 SUBROUTINE cp_cfm_diag_gen_dlaf(amatrix, bmatrix, eigenvectors, eigenvalues)
116
117 TYPE(cp_cfm_type), INTENT(IN) :: amatrix, bmatrix, eigenvectors
118 REAL(kind=dp), DIMENSION(:), INTENT(OUT) :: eigenvalues
119
120 CHARACTER(LEN=*), PARAMETER :: routinen = 'cp_cfm_diag_gen_dlaf'
121
122 INTEGER :: handle, n, nmo
123 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), TARGET :: eig
124
125 CALL timeset(routinen, handle)
126
127 n = amatrix%matrix_struct%nrow_global
128 ALLOCATE (eig(n))
129
130 CALL cp_cfm_diag_gen_dlaf_base(amatrix, bmatrix, eigenvectors, eig)
131
132 nmo = SIZE(eigenvalues, 1)
133 IF (nmo > n) THEN
134 eigenvalues(1:n) = eig(1:n)
135 ELSE
136 eigenvalues(1:nmo) = eig(1:nmo)
137 END IF
138
139 DEALLOCATE (eig)
140
141 CALL timestop(handle)
142
143 END SUBROUTINE cp_cfm_diag_gen_dlaf
144
145 !***************************************************************************************************
146!> \brief DLA-Future standard eigensolver for complex Hermitian matrices
147!> \param matrix ...
148!> \param eigenvectors ...
149!> \param eigenvalues ...
150!> \author Rocco Meli
151! **************************************************************************************************
152 SUBROUTINE cp_cfm_diag_dlaf_base(matrix, eigenvectors, eigenvalues)
153 TYPE(cp_cfm_type), INTENT(IN) :: matrix, eigenvectors
154 REAL(kind=dp), DIMENSION(:), INTENT(OUT), TARGET :: eigenvalues
155
156 CHARACTER(len=*), PARAMETER :: dlaf_name = 'pzheevd_dlaf', routinen = 'cp_cfm_diag_dlaf_base'
157 CHARACTER, PARAMETER :: uplo = 'L'
158
159 CHARACTER(LEN=100) :: message
160 COMPLEX(KIND=dp), DIMENSION(:, :), POINTER :: a, z
161 INTEGER :: dlaf_handle, handle, n
162 INTEGER, DIMENSION(9) :: desca, descz
163 INTEGER, TARGET :: info
164
165 CALL timeset(routinen, handle)
166
167#if defined(__DLAF)
168 ! DLAF needs the lower triangular part
169 ! Use eigenvectors matrix as workspace
170 CALL cp_cfm_uplo_to_full(matrix, eigenvectors)
171
172 n = matrix%matrix_struct%nrow_global
173
174 a => matrix%local_data
175 z => eigenvectors%local_data
176
177 desca(:) = matrix%matrix_struct%descriptor(:)
178 descz(:) = eigenvectors%matrix_struct%descriptor(:)
179
180 info = -1
181 CALL timeset(dlaf_name, dlaf_handle)
182
183 CALL dlaf_pzheevd(uplo, n, a, 1, 1, desca, eigenvalues, z, 1, 1, descz, info)
184
185 CALL timestop(dlaf_handle)
186
187 IF (info /= 0) THEN
188 WRITE (message, "(A,I0,A)") "ERROR in DLAF_PZHEEVD: Eigensolver failed (INFO = ", info, ")"
189 cpabort(trim(message))
190 END IF
191#else
192 mark_used(a)
193 mark_used(z)
194 mark_used(desca)
195 mark_used(descz)
196 mark_used(matrix)
197 mark_used(eigenvectors)
198 mark_used(eigenvalues)
199 mark_used(uplo)
200 mark_used(n)
201 mark_used(info)
202 mark_used(dlaf_handle)
203 mark_used(dlaf_name)
204 mark_used(message)
205 cpabort("CP2K compiled without DLAF library.")
206#endif
207
208 CALL timestop(handle)
209
210 END SUBROUTINE cp_cfm_diag_dlaf_base
211
212!***************************************************************************************************
213!> \brief DLA-Future generalized eigensolver for complex Hermitian matrices
214!> \param amatrix ...
215!> \param bmatrix ...
216!> \param eigenvectors ...
217!> \param eigenvalues ...
218!> \author Rocco Meli
219! **************************************************************************************************
220 SUBROUTINE cp_cfm_diag_gen_dlaf_base(amatrix, bmatrix, eigenvectors, eigenvalues)
221 TYPE(cp_cfm_type), INTENT(IN) :: amatrix, bmatrix, eigenvectors
222 REAL(kind=dp), DIMENSION(:), INTENT(OUT), TARGET :: eigenvalues
223
224 CHARACTER(len=*), PARAMETER :: dlaf_name = 'pzhegvd_dlaf', &
225 routinen = 'cp_cfm_diag_gen_dlaf_base'
226 CHARACTER, PARAMETER :: uplo = 'L'
227
228 CHARACTER(LEN=100) :: message
229 COMPLEX(KIND=dp), DIMENSION(:, :), POINTER :: a, b, z
230 INTEGER :: dlaf_handle, handle, n
231 INTEGER, DIMENSION(9) :: desca, descb, descz
232 INTEGER, TARGET :: info
233
234 CALL timeset(routinen, handle)
235
236#if defined(__DLAF)
237 ! DLAF needs the lower triangular part
238 ! Use eigenvectors matrix as workspace
239 CALL cp_cfm_uplo_to_full(amatrix, eigenvectors)
240 CALL cp_cfm_uplo_to_full(bmatrix, eigenvectors)
241
242 n = amatrix%matrix_struct%nrow_global
243
244 a => amatrix%local_data
245 b => bmatrix%local_data
246 z => eigenvectors%local_data
247
248 desca(:) = amatrix%matrix_struct%descriptor(:)
249 descb(:) = bmatrix%matrix_struct%descriptor(:)
250 descz(:) = eigenvectors%matrix_struct%descriptor(:)
251
252 info = -1
253 CALL timeset(dlaf_name, dlaf_handle)
254
255 CALL dlaf_pzhegvd(uplo, n, a, 1, 1, desca, b, 1, 1, descb, eigenvalues, z, 1, 1, descz, info)
256
257 CALL timestop(dlaf_handle)
258
259 IF (info /= 0) THEN
260 WRITE (message, "(A,I0,A)") "ERROR in DLAF_PZHEGVD: Eigensolver failed (INFO = ", info, ")"
261 cpabort(trim(message))
262 END IF
263#else
264 mark_used(a)
265 mark_used(b)
266 mark_used(z)
267 mark_used(desca)
268 mark_used(descb)
269 mark_used(descz)
270 mark_used(amatrix)
271 mark_used(bmatrix)
272 mark_used(eigenvectors)
273 mark_used(eigenvalues)
274 mark_used(uplo)
275 mark_used(n)
276 mark_used(info)
277 mark_used(dlaf_handle)
278 mark_used(dlaf_name)
279 mark_used(message)
280 cpabort("CP2K compiled without DLAF library.")
281#endif
282
283 CALL timestop(handle)
284
285 END SUBROUTINE cp_cfm_diag_gen_dlaf_base
286
287END MODULE cp_cfm_dlaf_api
Basic linear algebra operations for complex full matrices.
subroutine, public cp_cfm_uplo_to_full(matrix, workspace, uplo)
...
subroutine, public cp_cfm_diag_dlaf(matrix, eigenvectors, eigenvalues)
DLA-Future eigensolver for complex Hermitian matrices.
subroutine, public cp_cfm_pzpotrf_dlaf(uplo, n, a, ia, ja, desca, info)
Cholesky factorization using DLA-Future.
subroutine, public cp_cfm_diag_gen_dlaf(amatrix, bmatrix, eigenvectors, eigenvalues)
DLA-Future generalized eigensolver for complex Hermitian matrices.
Represents a complex full matrix distributed on many processors.
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Represent a complex full matrix.