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