(git:97501a3)
Loading...
Searching...
No Matches
parallel_gemm_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
8! **************************************************************************************************
9!> \brief basic linear algebra operations for full matrixes
10!> \par History
11!> 08.2002 splitted out of qs_blacs [fawzi]
12!> \author Fawzi Mohamed
13! **************************************************************************************************
15 USE iso_c_binding, ONLY: c_char,&
16 c_double,&
17 c_int,&
18 c_loc,&
19 c_ptr
21 USE cp_cfm_types, ONLY: cp_cfm_type
25 USE input_constants, ONLY: do_cosma,&
27 USE kinds, ONLY: dp
29#include "./base/base_uses.f90"
30
31 IMPLICIT NONE
32 PRIVATE
33
34 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'parallel_gemm_api'
35
36 PUBLIC :: parallel_gemm
37
38 INTERFACE parallel_gemm
39 MODULE PROCEDURE parallel_gemm_fm
40 MODULE PROCEDURE parallel_gemm_cfm
41 END INTERFACE parallel_gemm
42
43CONTAINS
44
45! **************************************************************************************************
46!> \brief ...
47!> \param transa ...
48!> \param transb ...
49!> \param m ...
50!> \param n ...
51!> \param k ...
52!> \param alpha ...
53!> \param matrix_a ...
54!> \param matrix_b ...
55!> \param beta ...
56!> \param matrix_c ...
57!> \param a_first_col ...
58!> \param a_first_row ...
59!> \param b_first_col ...
60!> \param b_first_row ...
61!> \param c_first_col ...
62!> \param c_first_row ...
63! **************************************************************************************************
64 SUBROUTINE parallel_gemm_fm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, &
65 matrix_c, a_first_col, a_first_row, b_first_col, b_first_row, &
66 c_first_col, c_first_row)
67 CHARACTER(LEN=1), INTENT(IN) :: transa, transb
68 INTEGER, INTENT(IN) :: m, n, k
69 REAL(KIND=dp), INTENT(IN) :: alpha
70 TYPE(cp_fm_type), INTENT(IN) :: matrix_a, matrix_b
71 REAL(kind=dp), INTENT(IN) :: beta
72 TYPE(cp_fm_type), INTENT(IN) :: matrix_c
73 INTEGER, INTENT(IN), OPTIONAL :: a_first_col, a_first_row, b_first_col, &
74 b_first_row, c_first_col, c_first_row
75
76 CHARACTER(len=*), PARAMETER :: routinen = 'parallel_gemm_fm'
77
78 INTEGER :: handle, my_multi
79
80 my_multi = cp_fm_get_mm_type()
81
82 SELECT CASE (my_multi)
83 CASE (do_scalapack)
84 CALL timeset(routinen//"_gemm", handle)
85 CALL cp_fm_gemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, &
86 a_first_col=a_first_col, &
87 a_first_row=a_first_row, &
88 b_first_col=b_first_col, &
89 b_first_row=b_first_row, &
90 c_first_col=c_first_col, &
91 c_first_row=c_first_row)
92 CASE (do_cosma)
93#if defined(__COSMA)
94 CALL timeset(routinen//"_cosma", handle)
96 CALL cosma_pdgemm(transa=transa, transb=transb, m=m, n=n, k=k, alpha=alpha, &
97 matrix_a=matrix_a, matrix_b=matrix_b, beta=beta, matrix_c=matrix_c, &
98 a_first_col=a_first_col, &
99 a_first_row=a_first_row, &
100 b_first_col=b_first_col, &
101 b_first_row=b_first_row, &
102 c_first_col=c_first_col, &
103 c_first_row=c_first_row)
104#else
105 cpabort("CP2K compiled without the COSMA library.")
106#endif
107 END SELECT
108 CALL timestop(handle)
109
110 END SUBROUTINE parallel_gemm_fm
111
112! **************************************************************************************************
113!> \brief ...
114!> \param transa ...
115!> \param transb ...
116!> \param m ...
117!> \param n ...
118!> \param k ...
119!> \param alpha ...
120!> \param matrix_a ...
121!> \param matrix_b ...
122!> \param beta ...
123!> \param matrix_c ...
124!> \param a_first_col ...
125!> \param a_first_row ...
126!> \param b_first_col ...
127!> \param b_first_row ...
128!> \param c_first_col ...
129!> \param c_first_row ...
130! **************************************************************************************************
131 SUBROUTINE parallel_gemm_cfm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, &
132 matrix_c, a_first_col, a_first_row, b_first_col, b_first_row, &
133 c_first_col, c_first_row)
134 CHARACTER(LEN=1), INTENT(IN) :: transa, transb
135 INTEGER, INTENT(IN) :: m, n, k
136 COMPLEX(KIND=dp), INTENT(IN) :: alpha
137 TYPE(cp_cfm_type), INTENT(IN) :: matrix_a, matrix_b
138 COMPLEX(KIND=dp), INTENT(IN) :: beta
139 TYPE(cp_cfm_type), INTENT(IN) :: matrix_c
140 INTEGER, INTENT(IN), OPTIONAL :: a_first_col, a_first_row, b_first_col, &
141 b_first_row, c_first_col, c_first_row
142
143 CHARACTER(len=*), PARAMETER :: routineN = 'parallel_gemm_cfm'
144
145 INTEGER :: handle, handle1, my_multi
146
147 CALL timeset(routinen, handle)
148
149 my_multi = cp_fm_get_mm_type()
150
151 SELECT CASE (my_multi)
152 CASE (do_scalapack)
153 CALL timeset(routinen//"_gemm", handle1)
154 CALL cp_cfm_gemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, &
155 a_first_col=a_first_col, &
156 a_first_row=a_first_row, &
157 b_first_col=b_first_col, &
158 b_first_row=b_first_row, &
159 c_first_col=c_first_col, &
160 c_first_row=c_first_row)
161 CALL timestop(handle1)
162 CASE (do_cosma)
163#if defined(__COSMA)
164 CALL timeset(routinen//"_cosma", handle1)
166 CALL cosma_pzgemm(transa=transa, transb=transb, m=m, n=n, k=k, alpha=alpha, &
167 matrix_a=matrix_a, matrix_b=matrix_b, beta=beta, matrix_c=matrix_c, &
168 a_first_col=a_first_col, &
169 a_first_row=a_first_row, &
170 b_first_col=b_first_col, &
171 b_first_row=b_first_row, &
172 c_first_col=c_first_col, &
173 c_first_row=c_first_row)
174 CALL timestop(handle1)
175#else
176 cpabort("CP2K compiled without the COSMA library.")
177#endif
178 END SELECT
179 CALL timestop(handle)
180
181 END SUBROUTINE parallel_gemm_cfm
182
183#if defined(__COSMA)
184! **************************************************************************************************
185!> \brief Fortran wrapper for cosma_pdgemm.
186!> \param transa ...
187!> \param transb ...
188!> \param m ...
189!> \param n ...
190!> \param k ...
191!> \param alpha ...
192!> \param matrix_a ...
193!> \param matrix_b ...
194!> \param beta ...
195!> \param matrix_c ...
196!> \param a_first_col ...
197!> \param a_first_row ...
198!> \param b_first_col ...
199!> \param b_first_row ...
200!> \param c_first_col ...
201!> \param c_first_row ...
202!> \author Ole Schuett
203! **************************************************************************************************
204 SUBROUTINE cosma_pdgemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, &
205 a_first_col, a_first_row, b_first_col, b_first_row, &
206 c_first_col, c_first_row)
207 CHARACTER(LEN=1), INTENT(IN) :: transa, transb
208 INTEGER, INTENT(IN) :: m, n, k
209 REAL(KIND=dp), INTENT(IN) :: alpha
210 TYPE(cp_fm_type), INTENT(IN) :: matrix_a, matrix_b
211 REAL(kind=dp), INTENT(IN) :: beta
212 TYPE(cp_fm_type), INTENT(IN) :: matrix_c
213 INTEGER, INTENT(IN), OPTIONAL :: a_first_col, a_first_row, b_first_col, &
214 b_first_row, c_first_col, c_first_row
215
216 INTEGER :: i_a, i_b, i_c, j_a, j_b, j_c
217 INTERFACE
218 SUBROUTINE cosma_pdgemm_c(transa, transb, m, n, k, alpha, a, ia, ja, desca, &
219 b, ib, jb, descb, beta, c, ic, jc, descc) &
220 BIND(C, name="cosma_pdgemm")
221 IMPORT :: c_ptr, c_int, c_double, c_char
222 CHARACTER(KIND=C_CHAR) :: transa
223 CHARACTER(KIND=C_CHAR) :: transb
224 INTEGER(KIND=C_INT) :: m
225 INTEGER(KIND=C_INT) :: n
226 INTEGER(KIND=C_INT) :: k
227 REAL(kind=c_double) :: alpha
228 TYPE(c_ptr), VALUE :: a
229 INTEGER(KIND=C_INT) :: ia
230 INTEGER(KIND=C_INT) :: ja
231 TYPE(c_ptr), VALUE :: desca
232 TYPE(c_ptr), VALUE :: b
233 INTEGER(KIND=C_INT) :: ib
234 INTEGER(KIND=C_INT) :: jb
235 TYPE(c_ptr), VALUE :: descb
236 REAL(KIND=c_double) :: beta
237 TYPE(c_ptr), VALUE :: c
238 INTEGER(KIND=C_INT) :: ic
239 INTEGER(KIND=C_INT) :: jc
240 TYPE(c_ptr), VALUE :: descc
241 END SUBROUTINE cosma_pdgemm_c
242 END INTERFACE
243
244 IF (PRESENT(a_first_row)) THEN
245 i_a = a_first_row
246 ELSE
247 i_a = 1
248 END IF
249 IF (PRESENT(a_first_col)) THEN
250 j_a = a_first_col
251 ELSE
252 j_a = 1
253 END IF
254 IF (PRESENT(b_first_row)) THEN
255 i_b = b_first_row
256 ELSE
257 i_b = 1
258 END IF
259 IF (PRESENT(b_first_col)) THEN
260 j_b = b_first_col
261 ELSE
262 j_b = 1
263 END IF
264 IF (PRESENT(c_first_row)) THEN
265 i_c = c_first_row
266 ELSE
267 i_c = 1
268 END IF
269 IF (PRESENT(c_first_col)) THEN
270 j_c = c_first_col
271 ELSE
272 j_c = 1
273 END IF
274
275 CALL cosma_pdgemm_c(transa=transa, transb=transb, m=m, n=n, k=k, &
276 alpha=alpha, &
277 a=c_loc(matrix_a%local_data(1, 1)), ia=i_a, ja=j_a, &
278 desca=c_loc(matrix_a%matrix_struct%descriptor(1)), &
279 b=c_loc(matrix_b%local_data(1, 1)), ib=i_b, jb=j_b, &
280 descb=c_loc(matrix_b%matrix_struct%descriptor(1)), &
281 beta=beta, &
282 c=c_loc(matrix_c%local_data(1, 1)), ic=i_c, jc=j_c, &
283 descc=c_loc(matrix_c%matrix_struct%descriptor(1)))
284
285 END SUBROUTINE cosma_pdgemm
286
287! **************************************************************************************************
288!> \brief Fortran wrapper for cosma_pdgemm.
289!> \param transa ...
290!> \param transb ...
291!> \param m ...
292!> \param n ...
293!> \param k ...
294!> \param alpha ...
295!> \param matrix_a ...
296!> \param matrix_b ...
297!> \param beta ...
298!> \param matrix_c ...
299!> \param a_first_col ...
300!> \param a_first_row ...
301!> \param b_first_col ...
302!> \param b_first_row ...
303!> \param c_first_col ...
304!> \param c_first_row ...
305!> \author Ole Schuett
306! **************************************************************************************************
307 SUBROUTINE cosma_pzgemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, &
308 a_first_col, a_first_row, b_first_col, b_first_row, &
309 c_first_col, c_first_row)
310 CHARACTER(LEN=1), INTENT(IN) :: transa, transb
311 INTEGER, INTENT(IN) :: m, n, k
312 COMPLEX(KIND=dp), INTENT(IN) :: alpha
313 TYPE(cp_cfm_type), INTENT(IN) :: matrix_a, matrix_b
314 COMPLEX(KIND=dp), INTENT(IN) :: beta
315 TYPE(cp_cfm_type), INTENT(IN) :: matrix_c
316 INTEGER, INTENT(IN), OPTIONAL :: a_first_col, a_first_row, b_first_col, &
317 b_first_row, c_first_col, c_first_row
318
319 INTEGER :: i_a, i_b, i_c, j_a, j_b, j_c
320 REAL(KIND=dp), DIMENSION(2), TARGET :: alpha_t, beta_t
321 INTERFACE
322 SUBROUTINE cosma_pzgemm_c(transa, transb, m, n, k, alpha, a, ia, ja, desca, &
323 b, ib, jb, descb, beta, c, ic, jc, descc) &
324 BIND(C, name="cosma_pzgemm")
325 IMPORT :: c_ptr, c_int, c_char
326 CHARACTER(KIND=C_CHAR) :: transa
327 CHARACTER(KIND=C_CHAR) :: transb
328 INTEGER(KIND=C_INT) :: m
329 INTEGER(KIND=C_INT) :: n
330 INTEGER(KIND=C_INT) :: k
331 TYPE(c_ptr), VALUE :: alpha
332 TYPE(c_ptr), VALUE :: a
333 INTEGER(KIND=C_INT) :: ia
334 INTEGER(KIND=C_INT) :: ja
335 TYPE(c_ptr), VALUE :: desca
336 TYPE(c_ptr), VALUE :: b
337 INTEGER(KIND=C_INT) :: ib
338 INTEGER(KIND=C_INT) :: jb
339 TYPE(c_ptr), VALUE :: descb
340 TYPE(c_ptr), VALUE :: beta
341 TYPE(c_ptr), VALUE :: c
342 INTEGER(KIND=C_INT) :: ic
343 INTEGER(KIND=C_INT) :: jc
344 TYPE(c_ptr), VALUE :: descc
345 END SUBROUTINE cosma_pzgemm_c
346 END INTERFACE
347
348 IF (PRESENT(a_first_row)) THEN
349 i_a = a_first_row
350 ELSE
351 i_a = 1
352 END IF
353 IF (PRESENT(a_first_col)) THEN
354 j_a = a_first_col
355 ELSE
356 j_a = 1
357 END IF
358 IF (PRESENT(b_first_row)) THEN
359 i_b = b_first_row
360 ELSE
361 i_b = 1
362 END IF
363 IF (PRESENT(b_first_col)) THEN
364 j_b = b_first_col
365 ELSE
366 j_b = 1
367 END IF
368 IF (PRESENT(c_first_row)) THEN
369 i_c = c_first_row
370 ELSE
371 i_c = 1
372 END IF
373 IF (PRESENT(c_first_col)) THEN
374 j_c = c_first_col
375 ELSE
376 j_c = 1
377 END IF
378
379 alpha_t(1) = real(alpha, kind=dp)
380 alpha_t(2) = real(aimag(alpha), kind=dp)
381 beta_t(1) = real(beta, kind=dp)
382 beta_t(2) = real(aimag(beta), kind=dp)
383
384 CALL cosma_pzgemm_c(transa=transa, transb=transb, m=m, n=n, k=k, &
385 alpha=c_loc(alpha_t), &
386 a=c_loc(matrix_a%local_data(1, 1)), ia=i_a, ja=j_a, &
387 desca=c_loc(matrix_a%matrix_struct%descriptor(1)), &
388 b=c_loc(matrix_b%local_data(1, 1)), ib=i_b, jb=j_b, &
389 descb=c_loc(matrix_b%matrix_struct%descriptor(1)), &
390 beta=c_loc(beta_t), &
391 c=c_loc(matrix_c%local_data(1, 1)), ic=i_c, jc=j_c, &
392 descc=c_loc(matrix_c%matrix_struct%descriptor(1)))
393
394 END SUBROUTINE cosma_pzgemm
395#endif
396
397END MODULE parallel_gemm_api
Basic linear algebra operations for complex full matrices.
subroutine, public cp_cfm_gemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, a_first_col, a_first_row, b_first_col, b_first_row, c_first_col, c_first_row)
Performs one of the matrix-matrix operations: matrix_c = alpha * op1( matrix_a ) * op2( matrix_b ) + ...
Represents a complex full matrix distributed on many processors.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_gemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, a_first_col, a_first_row, b_first_col, b_first_row, c_first_col, c_first_row)
computes matrix_c = beta * matrix_c + alpha * ( matrix_a ** transa ) * ( matrix_b ** transb )
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
integer function, public cp_fm_get_mm_type()
...
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_cosma
integer, parameter, public do_scalapack
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Fortran API for the offload package, which is written in C.
Definition offload_api.F:12
subroutine, public offload_activate_chosen_device()
Activates the device selected via offload_set_chosen_device()
basic linear algebra operations for full matrixes
Represent a complex full matrix.
represent a full matrix