Line data Source code
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 :
8 : MODULE cp_fm_dlaf_api
9 :
10 : USE cp_fm_basic_linalg, ONLY: cp_fm_uplo_to_full
11 : USE cp_fm_types, ONLY: cp_fm_type
12 : USE kinds, ONLY: sp, dp
13 : #include "../base/base_uses.f90"
14 :
15 : #if defined(__DLAF)
16 : USE cp_dlaf_utils_api, ONLY: cp_dlaf_create_grid
17 : USE dlaf_fortran, ONLY: dlaf_pdpotrf, &
18 : dlaf_pdsyevd, &
19 : dlaf_pdsygvd, &
20 : dlaf_pspotrf, &
21 : dlaf_pspotri, &
22 : dlaf_pdpotri
23 : #endif
24 :
25 : IMPLICIT NONE
26 :
27 : PRIVATE
28 :
29 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp_fm_dlaf_api'
30 :
31 : PUBLIC :: cp_pdpotrf_dlaf, cp_pspotrf_dlaf
32 : PUBLIC :: cp_pdpotri_dlaf, cp_pspotri_dlaf
33 : PUBLIC :: cp_fm_diag_dlaf, cp_fm_diag_gen_dlaf
34 :
35 : CONTAINS
36 :
37 : !***************************************************************************************************
38 : !> \brief Cholesky factorization using DLA-Future
39 : !> \param uplo ...
40 : !> \param n Matrix size
41 : !> \param a Local matrix
42 : !> \param ia Row index of first row (has to be 1)
43 : !> \param ja Col index of first column ()
44 : !> \param desca ScaLAPACK matrix descriptor
45 : !> \param info 0 if factorization completed normally
46 : !> \author Rocco Meli
47 : !> \author Mikael Simberg
48 : !> \author Mathieu Taillefumier
49 : ! **************************************************************************************************
50 0 : SUBROUTINE cp_pdpotrf_dlaf(uplo, n, a, ia, ja, desca, info)
51 : CHARACTER, INTENT(IN) :: uplo
52 : INTEGER, INTENT(IN) :: n
53 : REAL(KIND=dp), DIMENSION(:, :), TARGET :: a
54 : INTEGER, INTENT(IN) :: ia, ja
55 : INTEGER, DIMENSION(9) :: desca
56 : INTEGER, TARGET :: info
57 :
58 : CHARACTER(len=*), PARAMETER :: routineN = 'cp_pdpotrf_dlaf'
59 :
60 : INTEGER :: handle
61 :
62 0 : CALL timeset(routineN, handle)
63 : #if defined(__DLAF)
64 : CALL dlaf_pdpotrf(uplo, n, a, ia, ja, desca, info)
65 : #else
66 : MARK_USED(uplo)
67 : MARK_USED(n)
68 : MARK_USED(a)
69 : MARK_USED(ia)
70 : MARK_USED(ja)
71 : MARK_USED(desca)
72 : MARK_USED(info)
73 0 : CPABORT("CP2K compiled without the DLA-Future library.")
74 : #endif
75 0 : CALL timestop(handle)
76 0 : END SUBROUTINE cp_pdpotrf_dlaf
77 :
78 : !***************************************************************************************************
79 : !> \brief Cholesky factorization using DLA-Future
80 : !> \param uplo ...
81 : !> \param n Matrix size
82 : !> \param a Local matrix
83 : !> \param ia Row index of first row (has to be 1)
84 : !> \param ja Col index of first column ()
85 : !> \param desca ScaLAPACK matrix descriptor
86 : !> \param info 0 if factorization completed normally
87 : !> \author Rocco Meli
88 : !> \author Mikael Simberg
89 : !> \author Mathieu Taillefumier
90 : ! **************************************************************************************************
91 0 : SUBROUTINE cp_pspotrf_dlaf(uplo, n, a, ia, ja, desca, info)
92 : CHARACTER, INTENT(IN) :: uplo
93 : INTEGER, INTENT(IN) :: n
94 : REAL, DIMENSION(:, :), TARGET :: a
95 : INTEGER, INTENT(IN) :: ia, ja
96 : INTEGER, DIMENSION(9) :: desca
97 : INTEGER, TARGET :: info
98 :
99 : CHARACTER(len=*), PARAMETER :: routineN = 'cp_pspotrf_dlaf'
100 :
101 : INTEGER :: handle
102 :
103 0 : CALL timeset(routineN, handle)
104 : #if defined(__DLAF)
105 : CALL dlaf_pspotrf(uplo, n, a, ia, ja, desca, info)
106 : #else
107 : MARK_USED(uplo)
108 : MARK_USED(n)
109 : MARK_USED(a)
110 : MARK_USED(ia)
111 : MARK_USED(ja)
112 : MARK_USED(desca)
113 : MARK_USED(info)
114 0 : CPABORT("CP2K compiled without the DLA-Future library.")
115 : #endif
116 0 : CALL timestop(handle)
117 0 : END SUBROUTINE cp_pspotrf_dlaf
118 :
119 : !***************************************************************************************************
120 : !> \brief Inverse from Cholesky factorization using DLA-Future
121 : !> \param uplo ...
122 : !> \param n Matrix size
123 : !> \param a Local matrix
124 : !> \param ia Row index of first row (has to be 1)
125 : !> \param ja Col index of first column ()
126 : !> \param desca ScaLAPACK matrix descriptor
127 : !> \param info 0 if factorization completed normally
128 : !> \author Rocco Meli
129 : ! **************************************************************************************************
130 0 : SUBROUTINE cp_pdpotri_dlaf(uplo, n, a, ia, ja, desca, info)
131 : CHARACTER, INTENT(IN) :: uplo
132 : INTEGER, INTENT(IN) :: n
133 : REAL(KIND=dp), DIMENSION(:, :), TARGET :: a
134 : INTEGER, INTENT(IN) :: ia, ja
135 : INTEGER, DIMENSION(9) :: desca
136 : INTEGER, TARGET :: info
137 :
138 : CHARACTER(len=*), PARAMETER :: routineN = 'cp_pdpotri_dlaf'
139 :
140 : INTEGER :: handle
141 :
142 0 : CALL timeset(routineN, handle)
143 : #if defined(__DLAF)
144 : CALL dlaf_pdpotri(uplo, n, a, ia, ja, desca, info)
145 : #else
146 : MARK_USED(uplo)
147 : MARK_USED(n)
148 : MARK_USED(a)
149 : MARK_USED(ia)
150 : MARK_USED(ja)
151 : MARK_USED(desca)
152 : MARK_USED(info)
153 0 : CPABORT("CP2K compiled without the DLA-Future library.")
154 : #endif
155 0 : CALL timestop(handle)
156 0 : END SUBROUTINE cp_pdpotri_dlaf
157 :
158 : !***************************************************************************************************
159 : !> \brief Inverse from Cholesky factorization using DLA-Future
160 : !> \param uplo ...
161 : !> \param n Matrix size
162 : !> \param a Local matrix
163 : !> \param ia Row index of first row (has to be 1)
164 : !> \param ja Col index of first column ()
165 : !> \param desca ScaLAPACK matrix descriptor
166 : !> \param info 0 if factorization completed normally
167 : !> \author Rocco Meli
168 : ! **************************************************************************************************
169 0 : SUBROUTINE cp_pspotri_dlaf(uplo, n, a, ia, ja, desca, info)
170 : CHARACTER, INTENT(IN) :: uplo
171 : INTEGER, INTENT(IN) :: n
172 : REAL, DIMENSION(:, :), TARGET :: a
173 : INTEGER, INTENT(IN) :: ia, ja
174 : INTEGER, DIMENSION(9) :: desca
175 : INTEGER, TARGET :: info
176 :
177 : CHARACTER(len=*), PARAMETER :: routineN = 'cp_pspotri_dlaf'
178 :
179 : INTEGER :: handle
180 :
181 0 : CALL timeset(routineN, handle)
182 :
183 : #if defined(__DLAF)
184 : CALL dlaf_pspotri(uplo, n, a, ia, ja, desca, info)
185 : #else
186 : MARK_USED(uplo)
187 : MARK_USED(n)
188 : MARK_USED(a)
189 : MARK_USED(ia)
190 : MARK_USED(ja)
191 : MARK_USED(desca)
192 : MARK_USED(info)
193 0 : CPABORT("CP2K compiled without the DLA-Future library.")
194 : #endif
195 0 : CALL timestop(handle)
196 0 : END SUBROUTINE cp_pspotri_dlaf
197 :
198 : ! **************************************************************************************************
199 : !> \brief ...
200 : !> \param matrix ...
201 : !> \param eigenvectors ...
202 : !> \param eigenvalues ...
203 : ! **************************************************************************************************
204 0 : SUBROUTINE cp_fm_diag_dlaf(matrix, eigenvectors, eigenvalues)
205 :
206 : TYPE(cp_fm_type), INTENT(IN) :: matrix, eigenvectors
207 : REAL(KIND=dp), DIMENSION(:), INTENT(OUT) :: eigenvalues
208 :
209 : CHARACTER(LEN=*), PARAMETER :: routineN = 'cp_fm_diag_dlaf'
210 :
211 : INTEGER :: handle, n, nmo
212 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), TARGET :: eig
213 :
214 0 : CALL timeset(routineN, handle)
215 :
216 0 : n = matrix%matrix_struct%nrow_global
217 0 : ALLOCATE (eig(n))
218 :
219 0 : CALL cp_fm_diag_dlaf_base(matrix, eigenvectors, eig)
220 :
221 0 : nmo = SIZE(eigenvalues, 1)
222 0 : IF (nmo > n) THEN
223 0 : eigenvalues(1:n) = eig(1:n)
224 : ELSE
225 0 : eigenvalues(1:nmo) = eig(1:nmo)
226 : END IF
227 :
228 0 : DEALLOCATE (eig)
229 :
230 0 : CALL timestop(handle)
231 :
232 0 : END SUBROUTINE cp_fm_diag_dlaf
233 :
234 : !***************************************************************************************************
235 : !> \brief DLA-Future eigensolver
236 : !> \param matrix ...
237 : !> \param eigenvectors ...
238 : !> \param eigenvalues ...
239 : !> \author Rocco Meli
240 : ! **************************************************************************************************
241 0 : SUBROUTINE cp_fm_diag_dlaf_base(matrix, eigenvectors, eigenvalues)
242 : TYPE(cp_fm_type), INTENT(IN) :: matrix, eigenvectors
243 : REAL(kind=dp), DIMENSION(:), INTENT(OUT), TARGET :: eigenvalues
244 :
245 : CHARACTER(len=*), PARAMETER :: dlaf_name = 'pdsyevd_dlaf', routineN = 'cp_fm_diag_dlaf_base'
246 : CHARACTER, PARAMETER :: uplo = 'L'
247 :
248 : CHARACTER(LEN=100) :: message
249 : INTEGER :: blacs_context, dlaf_handle, handle, n
250 : INTEGER, DIMENSION(9) :: desca, descz
251 : INTEGER, TARGET :: info
252 0 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: a, z
253 :
254 0 : CALL timeset(routineN, handle)
255 :
256 : #if defined(__DLAF)
257 : ! DLAF needs the lower triangular part
258 : ! Use eigenvectors matrix as workspace
259 : CALL cp_fm_uplo_to_full(matrix, eigenvectors)
260 :
261 : ! Create DLAF grid from BLACS context (if already present, does nothing)
262 : blacs_context = matrix%matrix_struct%context%get_handle()
263 : CALL cp_dlaf_create_grid(blacs_context)
264 :
265 : n = matrix%matrix_struct%nrow_global
266 :
267 : a => matrix%local_data
268 : z => eigenvectors%local_data
269 :
270 : desca(:) = matrix%matrix_struct%descriptor(:)
271 : descz(:) = eigenvectors%matrix_struct%descriptor(:)
272 :
273 : info = -1
274 : CALL timeset(dlaf_name, dlaf_handle)
275 : CALL dlaf_pdsyevd(uplo, n, a, 1, 1, desca, eigenvalues, z, 1, 1, descz, info)
276 : CALL timestop(dlaf_handle)
277 :
278 : IF (info /= 0) THEN
279 : WRITE (message, "(A,I0,A)") "ERROR in DLAF_PDSYEVD: Eigensolver failed (INFO = ", info, ")"
280 : CPABORT(TRIM(message))
281 : END IF
282 : #else
283 : MARK_USED(a)
284 : MARK_USED(z)
285 : MARK_USED(desca)
286 : MARK_USED(descz)
287 : MARK_USED(matrix)
288 : MARK_USED(eigenvectors)
289 : MARK_USED(eigenvalues)
290 : MARK_USED(uplo)
291 : MARK_USED(n)
292 : MARK_USED(info)
293 : MARK_USED(dlaf_handle)
294 : MARK_USED(dlaf_name)
295 : MARK_USED(message)
296 : MARK_USED(blacs_context)
297 0 : CPABORT("CP2K compiled without DLA-Future-Fortran library.")
298 : #endif
299 :
300 0 : CALL timestop(handle)
301 :
302 0 : END SUBROUTINE cp_fm_diag_dlaf_base
303 :
304 : ! **************************************************************************************************
305 : !> \brief ...
306 : !> \param a_matrix ...
307 : !> \param b_matrix ...
308 : !> \param eigenvectors ...
309 : !> \param eigenvalues ...
310 : !> \author Rocco Meli
311 : ! **************************************************************************************************
312 0 : SUBROUTINE cp_fm_diag_gen_dlaf(a_matrix, b_matrix, eigenvectors, eigenvalues)
313 :
314 : TYPE(cp_fm_type), INTENT(IN) :: a_matrix, b_matrix, eigenvectors
315 : REAL(KIND=dp), DIMENSION(:), INTENT(OUT) :: eigenvalues
316 :
317 : CHARACTER(LEN=*), PARAMETER :: routineN = 'cp_fm_diag_gen_dlaf'
318 :
319 : INTEGER :: handle, n, nmo
320 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), TARGET :: eig
321 :
322 0 : CALL timeset(routineN, handle)
323 :
324 0 : n = a_matrix%matrix_struct%nrow_global
325 0 : ALLOCATE (eig(n))
326 :
327 0 : CALL cp_fm_diag_gen_dlaf_base(a_matrix, b_matrix, eigenvectors, eig)
328 :
329 0 : nmo = SIZE(eigenvalues, 1)
330 0 : IF (nmo > n) THEN
331 0 : eigenvalues(1:n) = eig(1:n)
332 : ELSE
333 0 : eigenvalues(1:nmo) = eig(1:nmo)
334 : END IF
335 :
336 0 : DEALLOCATE (eig)
337 :
338 0 : CALL timestop(handle)
339 :
340 0 : END SUBROUTINE cp_fm_diag_gen_dlaf
341 :
342 : !***************************************************************************************************
343 : !> \brief DLA-Future generalized eigensolver
344 : !> \param a_matrix ...
345 : !> \param b_matrix ...
346 : !> \param eigenvectors ...
347 : !> \param eigenvalues ...
348 : !> \author Rocco Meli
349 : ! **************************************************************************************************
350 0 : SUBROUTINE cp_fm_diag_gen_dlaf_base(a_matrix, b_matrix, eigenvectors, eigenvalues)
351 : TYPE(cp_fm_type), INTENT(IN) :: a_matrix, b_matrix, eigenvectors
352 : REAL(kind=dp), DIMENSION(:), INTENT(OUT), TARGET :: eigenvalues
353 :
354 : CHARACTER(len=*), PARAMETER :: dlaf_name = 'pdsyevd_dlaf', &
355 : routineN = 'cp_fm_diag_gen_dlaf_base'
356 : CHARACTER, PARAMETER :: uplo = 'L'
357 :
358 : CHARACTER(LEN=100) :: message
359 : INTEGER :: blacs_context, dlaf_handle, handle, n
360 : INTEGER, DIMENSION(9) :: desca, descb, descz
361 : INTEGER, TARGET :: info
362 0 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: a, b, z
363 :
364 0 : CALL timeset(routineN, handle)
365 :
366 : #if defined(__DLAF)
367 : ! DLAF needs the lower triangular part
368 : ! Use eigenvectors matrix as workspace
369 : CALL cp_fm_uplo_to_full(a_matrix, eigenvectors)
370 : CALL cp_fm_uplo_to_full(b_matrix, eigenvectors)
371 :
372 : ! Create DLAF grid from BLACS context; if already present, does nothing
373 : blacs_context = a_matrix%matrix_struct%context%get_handle()
374 : CALL cp_dlaf_create_grid(blacs_context)
375 :
376 : n = a_matrix%matrix_struct%nrow_global
377 :
378 : a => a_matrix%local_data
379 : b => b_matrix%local_data
380 : z => eigenvectors%local_data
381 :
382 : desca(:) = a_matrix%matrix_struct%descriptor(:)
383 : descb(:) = b_matrix%matrix_struct%descriptor(:)
384 : descz(:) = eigenvectors%matrix_struct%descriptor(:)
385 :
386 : info = -1
387 : CALL timeset(dlaf_name, dlaf_handle)
388 : CALL dlaf_pdsygvd(uplo, n, a, 1, 1, desca, b, 1, 1, descb, eigenvalues, z, 1, 1, descz, info)
389 : CALL timestop(dlaf_handle)
390 :
391 : IF (info /= 0) THEN
392 : WRITE (message, "(A,I0,A)") "ERROR in DLAF_PDSYGVD: Generalized Eigensolver failed (INFO = ", info, ")"
393 : CPABORT(TRIM(message))
394 : END IF
395 : #else
396 : MARK_USED(a)
397 : MARK_USED(b)
398 : MARK_USED(z)
399 : MARK_USED(desca)
400 : MARK_USED(descb)
401 : MARK_USED(descz)
402 : MARK_USED(a_matrix)
403 : MARK_USED(b_matrix)
404 : MARK_USED(eigenvectors)
405 : MARK_USED(eigenvalues)
406 : MARK_USED(uplo)
407 : MARK_USED(n)
408 : MARK_USED(info)
409 : MARK_USED(blacs_context)
410 : MARK_USED(dlaf_handle)
411 : MARK_USED(dlaf_name)
412 : MARK_USED(message)
413 0 : CPABORT("CP2K compiled without DLA-Future-Fortran library.")
414 : #endif
415 :
416 0 : CALL timestop(handle)
417 :
418 0 : END SUBROUTINE cp_fm_diag_gen_dlaf_base
419 :
420 : END MODULE cp_fm_dlaf_api
|