Line data Source code
1 : ! This file is part of ELPA.
2 : !
3 : ! The ELPA library was originally created by the ELPA consortium,
4 : ! consisting of the following organizations:
5 : !
6 : ! - Max Planck Computing and Data Facility (MPCDF), formerly known as
7 : ! Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
8 : ! - Bergische Universität Wuppertal, Lehrstuhl für angewandte
9 : ! Informatik,
10 : ! - Technische Universität München, Lehrstuhl für Informatik mit
11 : ! Schwerpunkt Wissenschaftliches Rechnen ,
12 : ! - Fritz-Haber-Institut, Berlin, Abt. Theorie,
13 : ! - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
14 : ! Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
15 : ! and
16 : ! - IBM Deutschland GmbH
17 : !
18 : !
19 : ! More information can be found here:
20 : ! http://elpa.mpcdf.mpg.de/
21 : !
22 : ! ELPA is free software: you can redistribute it and/or modify
23 : ! it under the terms of the version 3 of the license of the
24 : ! GNU Lesser General Public License as published by the Free
25 : ! Software Foundation.
26 : !
27 : ! ELPA is distributed in the hope that it will be useful,
28 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : ! GNU Lesser General Public License for more details.
31 : !
32 : ! You should have received a copy of the GNU Lesser General Public License
33 : ! along with ELPA. If not, see <http://www.gnu.org/licenses/>
34 : !
35 : ! ELPA reflects a substantial effort on the part of the original
36 : ! ELPA consortium, and we ask you to respect the spirit of the
37 : ! license that we chose: i.e., please contribute any changes you
38 : ! may have back to the original ELPA library distribution, and keep
39 : ! any derivatives of ELPA under the same license that we chose for
40 : ! the original distribution, the GNU Lesser General Public License.
41 : !
42 : ! Author: Andreas Marek, MCPDF
43 : #include "config-f90.h"
44 : !lc> #include <complex.h>
45 :
46 : !lc> /*! \brief C interface to create ELPA communicators
47 : !lc> *
48 : !lc> * \param mpi_comm_word MPI global communicator (in)
49 : !lc> * \param my_prow Row coordinate of the calling process in the process grid (in)
50 : !lc> * \param my_pcol Column coordinate of the calling process in the process grid (in)
51 : !lc> * \param mpi_comm_rows Communicator for communicating within rows of processes (out)
52 : !lc> * \result int integer error value of mpi_comm_split function
53 : !lc> */
54 : !lc> int elpa_get_communicators(int mpi_comm_world, int my_prow, int my_pcol, int *mpi_comm_rows, int *mpi_comm_cols);
55 2688 : function elpa_get_communicators_wrapper_c(mpi_comm_world, my_prow, my_pcol, &
56 : mpi_comm_rows, mpi_comm_cols) &
57 : result(mpierr) bind(C,name="elpa_get_communicators")
58 : use, intrinsic :: iso_c_binding
59 : use elpa1, only : elpa_get_communicators
60 :
61 : implicit none
62 : integer(kind=c_int) :: mpierr
63 : integer(kind=c_int), value :: mpi_comm_world, my_prow, my_pcol
64 : integer(kind=c_int) :: mpi_comm_rows, mpi_comm_cols
65 :
66 : mpierr = elpa_get_communicators(mpi_comm_world, my_prow, my_pcol, &
67 2688 : mpi_comm_rows, mpi_comm_cols)
68 :
69 2688 : end function
70 :
71 :
72 : !lc> /*! \brief C interface to solve the double-precision real eigenvalue problem with 1-stage solver
73 : !lc> *
74 : !lc> * \param na Order of matrix a
75 : !lc> * \param nev Number of eigenvalues needed.
76 : !lc> * The smallest nev eigenvalues/eigenvectors are calculated.
77 : !lc> * \param a Distributed matrix for which eigenvalues are to be computed.
78 : !lc> * Distribution is like in Scalapack.
79 : !lc> * The full matrix must be set (not only one half like in scalapack).
80 : !lc> * \param lda Leading dimension of a
81 : !lc> * \param ev(na) On output: eigenvalues of a, every processor gets the complete set
82 : !lc> * \param q On output: Eigenvectors of a
83 : !lc> * Distribution is like in Scalapack.
84 : !lc> * Must be always dimensioned to the full size (corresponding to (na,na))
85 : !lc> * even if only a part of the eigenvalues is needed.
86 : !lc> * \param ldq Leading dimension of q
87 : !lc> * \param nblk blocksize of cyclic distribution, must be the same in both directions!
88 : !lc> * \param matrixCols distributed number of matrix columns
89 : !lc> * \param mpi_comm_rows MPI-Communicator for rows
90 : !lc> * \param mpi_comm_cols MPI-Communicator for columns
91 : !lc> * \param useGPU use GPU (1=yes, 0=No)
92 : !lc> *
93 : !lc> * \result int: 1 if error occured, otherwise 0
94 : !lc>*/
95 : #define REALCASE 1
96 : #define DOUBLE_PRECISION 1
97 : #include "../../general/precision_macros.h"
98 : !lc> int elpa_solve_evp_real_1stage_double_precision(int na, int nev, double *a, int lda, double *ev, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
99 :
100 : #include "./elpa1_c_interface_template.F90"
101 : #undef REALCASE
102 : #undef DOUBLE_PRECISION
103 :
104 : #ifdef WANT_SINGLE_PRECISION_REAL
105 :
106 : !lc> /*! \brief C interface to solve the single-precision real eigenvalue problem with 1-stage solver
107 : !lc> *
108 : !lc> * \param na Order of matrix a
109 : !lc> * \param nev Number of eigenvalues needed.
110 : !lc> * The smallest nev eigenvalues/eigenvectors are calculated.
111 : !lc> * \param a Distributed matrix for which eigenvalues are to be computed.
112 : !lc> * Distribution is like in Scalapack.
113 : !lc> * The full matrix must be set (not only one half like in scalapack).
114 : !lc> * \param lda Leading dimension of a
115 : !lc> * \param ev(na) On output: eigenvalues of a, every processor gets the complete set
116 : !lc> * \param q On output: Eigenvectors of a
117 : !lc> * Distribution is like in Scalapack.
118 : !lc> * Must be always dimensioned to the full size (corresponding to (na,na))
119 : !lc> * even if only a part of the eigenvalues is needed.
120 : !lc> * \param ldq Leading dimension of q
121 : !lc> * \param nblk blocksize of cyclic distribution, must be the same in both directions!
122 : !lc> * \param matrixCols distributed number of matrix columns
123 : !lc> * \param mpi_comm_rows MPI-Communicator for rows
124 : !lc> * \param mpi_comm_cols MPI-Communicator for columns
125 : !lc> * \param useGPU use GPU (1=yes, 0=No)
126 : !lc> *
127 : !lc> * \result int: 1 if error occured, otherwise 0
128 : !lc>*/
129 : #define REALCASE 1
130 : #undef DOUBLE_PRECISION
131 : #define SINGLE_PRECISION 1
132 : #include "../../general/precision_macros.h"
133 :
134 : !lc> int elpa_solve_evp_real_1stage_single_precision(int na, int nev, float *a, int lda, float *ev, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
135 :
136 : #include "./elpa1_c_interface_template.F90"
137 : #undef SINGLE_PRECISION
138 : #undef REALCASE
139 : #endif /* WANT_SINGLE_PRECISION_REAL */
140 :
141 : !lc> /*! \brief C interface to solve the double-precision complex eigenvalue problem with 1-stage solver
142 : !lc> *
143 : !lc> * \param na Order of matrix a
144 : !lc> * \param nev Number of eigenvalues needed.
145 : !lc> * The smallest nev eigenvalues/eigenvectors are calculated.
146 : !lc> * \param a Distributed matrix for which eigenvalues are to be computed.
147 : !lc> * Distribution is like in Scalapack.
148 : !lc> * The full matrix must be set (not only one half like in scalapack).
149 : !lc> * \param lda Leading dimension of a
150 : !lc> * \param ev(na) On output: eigenvalues of a, every processor gets the complete set
151 : !lc> * \param q On output: Eigenvectors of a
152 : !lc> * Distribution is like in Scalapack.
153 : !lc> * Must be always dimensioned to the full size (corresponding to (na,na))
154 : !lc> * even if only a part of the eigenvalues is needed.
155 : !lc> * \param ldq Leading dimension of q
156 : !lc> * \param nblk blocksize of cyclic distribution, must be the same in both directions!
157 : !lc> * \param matrixCols distributed number of matrix columns
158 : !lc> * \param mpi_comm_rows MPI-Communicator for rows
159 : !lc> * \param mpi_comm_cols MPI-Communicator for columns
160 : !lc> * \param useGPU use GPU (1=yes, 0=No)
161 : !lc> *
162 : !lc> * \result int: 1 if error occured, otherwise 0
163 : !lc> */
164 :
165 : #define COMPLEXCASE 1
166 : #define DOUBLE_PRECISION 1
167 : #include "../../general/precision_macros.h"
168 :
169 : !lc> int elpa_solve_evp_complex_1stage_double_precision(int na, int nev, double complex *a, int lda, double *ev, double complex *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
170 :
171 : #include "./elpa1_c_interface_template.F90"
172 : #undef COMPLEXCASE
173 : #undef DOUBLE_PRECISION
174 :
175 : #ifdef WANT_SINGLE_PRECISION_COMPLEX
176 :
177 : !lc> /*! \brief C interface to solve the single-precision complex eigenvalue problem with 1-stage solver
178 : !lc> *
179 : !lc> * \param na Order of matrix a
180 : !lc> * \param nev Number of eigenvalues needed.
181 : !lc> * The smallest nev eigenvalues/eigenvectors are calculated.
182 : !lc> * \param a Distributed matrix for which eigenvalues are to be computed.
183 : !lc> * Distribution is like in Scalapack.
184 : !lc> * The full matrix must be set (not only one half like in scalapack).
185 : !lc> * \param lda Leading dimension of a
186 : !lc> * \param ev(na) On output: eigenvalues of a, every processor gets the complete set
187 : !lc> * \param q On output: Eigenvectors of a
188 : !lc> * Distribution is like in Scalapack.
189 : !lc> * Must be always dimensioned to the full size (corresponding to (na,na))
190 : !lc> * even if only a part of the eigenvalues is needed.
191 : !lc> * \param ldq Leading dimension of q
192 : !lc> * \param nblk blocksize of cyclic distribution, must be the same in both directions!
193 : !lc> * \param matrixCols distributed number of matrix columns
194 : !lc> * \param mpi_comm_rows MPI-Communicator for rows
195 : !lc> * \param mpi_comm_cols MPI-Communicator for columns
196 : !lc> * \param useGPU use GPU (1=yes, 0=No)
197 : !lc> *
198 : !lc> * \result int: 1 if error occured, otherwise 0
199 : !lc> */
200 : #define COMPLEXCASE 1
201 : #undef DOUBLE_PRECISION
202 : #define SINGLE_PRECISION
203 : #include "../../general/precision_macros.h"
204 :
205 : !lc> int elpa_solve_evp_complex_1stage_single_precision(int na, int nev, complex float *a, int lda, float *ev, complex float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int mpi_comm_all, int useGPU);
206 :
207 : #include "./elpa1_c_interface_template.F90"
208 :
209 : #undef SINGLE_PRECISION
210 : #undef COMPLEXCASE
211 : #endif /* WANT_SINGLE_PRECISION_COMPLEX */
212 :
213 : !lc> /*
214 : !lc> \brief C interface to solve double-precision tridiagonal eigensystem with divide and conquer method
215 : !lc> \details
216 : !lc>
217 : !lc> *\param na Matrix dimension
218 : !lc> *\param nev number of eigenvalues/vectors to be computed
219 : !lc> *\param d array d(na) on input diagonal elements of tridiagonal matrix, on
220 : !lc> * output the eigenvalues in ascending order
221 : !lc> *\param e array e(na) on input subdiagonal elements of matrix, on exit destroyed
222 : !lc> *\param q on exit : matrix q(ldq,matrixCols) contains the eigenvectors
223 : !lc> *\param ldq leading dimension of matrix q
224 : !lc> *\param nblk blocksize of cyclic distribution, must be the same in both directions!
225 : !lc> *\param matrixCols columns of matrix q
226 : !lc> *\param mpi_comm_rows MPI communicator for rows
227 : !lc> *\param mpi_comm_cols MPI communicator for columns
228 : !lc> *\param wantDebug give more debug information if 1, else 0
229 : !lc> *\result success int 1 on success, else 0
230 : !lc> */
231 : !lc> int elpa_solve_tridi_double(int na, int nev, double *d, double *e, double *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
232 : #define REALCASE 1
233 : #define DOUBLE_PRECISION 1
234 : #include "../../general/precision_macros.h"
235 : #include "./elpa_solve_tridi_c_interface_template.F90"
236 : #undef DOUBLE_PRECISION
237 : #undef REALCASE
238 :
239 : #ifdef WANT_SINGLE_PRECISION_REAL
240 :
241 : !lc> /*
242 : !lc> \brief C interface to solve single-precision tridiagonal eigensystem with divide and conquer method
243 : !lc> \details
244 : !lc>
245 : !lc> \param na Matrix dimension
246 : !lc> \param nev number of eigenvalues/vectors to be computed
247 : !lc> \param d array d(na) on input diagonal elements of tridiagonal matrix, on
248 : !lc> output the eigenvalues in ascending order
249 : !lc> \param e array e(na) on input subdiagonal elements of matrix, on exit destroyed
250 : !lc> \param q on exit : matrix q(ldq,matrixCols) contains the eigenvectors
251 : !lc> \param ldq leading dimension of matrix q
252 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
253 : !lc> \param matrixCols columns of matrix q
254 : !lc> \param mpi_comm_rows MPI communicator for rows
255 : !lc> \param mpi_comm_cols MPI communicator for columns
256 : !lc> \param wantDebug give more debug information if 1, else 0
257 : !lc> \result success int 1 on success, else 0
258 : !lc> */
259 : !lc> int elpa_solve_tridi_single(int na, int nev, float *d, float *e, float *q, int ldq, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
260 : #define REALCASE 1
261 : #define SINGLE_PRECISION 1
262 : #include "../../general/precision_macros.h"
263 : #include "./elpa_solve_tridi_c_interface_template.F90"
264 : #undef SINGLE_PRECISION
265 : #undef REALCASE
266 :
267 : #endif /* WANT_SINGLE_PRECISION_REAL */
268 :
269 : !lc> /*
270 : !lc> \brief C interface for elpa_mult_at_b_real_double: Performs C : = A**T * B for double-precision matrices
271 : !lc> where A is a square matrix (na,na) which is optionally upper or lower triangular
272 : !lc> B is a (na,ncb) matrix
273 : !lc> C is a (na,ncb) matrix where optionally only the upper or lower
274 : !lc> triangle may be computed
275 : !lc> \details
276 : !lc> \param uplo_a 'U' if A is upper triangular
277 : !lc> 'L' if A is lower triangular
278 : !lc> anything else if A is a full matrix
279 : !lc> Please note: This pertains to the original A (as set in the calling program)
280 : !lc> whereas the transpose of A is used for calculations
281 : !lc> If uplo_a is 'U' or 'L', the other triangle is not used at all,
282 : !lc> i.e. it may contain arbitrary numbers
283 : !lc> \param uplo_c 'U' if only the upper diagonal part of C is needed
284 : !lc> 'L' if only the upper diagonal part of C is needed
285 : !lc> anything else if the full matrix C is needed
286 : !lc> Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
287 : !lc> written to a certain extent, i.e. one shouldn't rely on the content there!
288 : !lc> \param na Number of rows/columns of A, number of rows of B and C
289 : !lc> \param ncb Number of columns of B and C
290 : !lc> \param a matrix a
291 : !lc> \param lda leading dimension of matrix a
292 : !lc> \param ldaCols columns of matrix a
293 : !lc> \param b matrix b
294 : !lc> \param ldb leading dimension of matrix b
295 : !lc> \param ldbCols columns of matrix b
296 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
297 : !lc> \param mpi_comm_rows MPI communicator for rows
298 : !lc> \param mpi_comm_cols MPI communicator for columns
299 : !lc> \param c matrix c
300 : !lc> \param ldc leading dimension of matrix c
301 : !lc> \param ldcCols columns of matrix c
302 : !lc> \result success int report success (1) or failure (0)
303 : !lc> */
304 :
305 : !lc> int elpa_mult_at_b_real_double(char uplo_a, char uplo_c, int na, int ncb, double *a, int lda, int ldaCols, double *b, int ldb, int ldbCols, int nlbk, int mpi_comm_rows, int mpi_comm_cols, double *c, int ldc, int ldcCols);
306 :
307 : #define REALCASE 1
308 : #define DOUBLE_PRECISION 1
309 : #include "../../general/precision_macros.h"
310 : #include "./elpa_mult_at_b_c_interface_template.F90"
311 : #undef DOUBLE_PRECISION
312 : #undef REALCASE
313 :
314 : #ifdef WANT_SINGLE_PRECISION_REAL
315 : !lc> /*
316 : !lc> \brief C interface for elpa_mult_at_b_real_single: Performs C : = A**T * B for single-precision matrices
317 : !lc> where A is a square matrix (na,na) which is optionally upper or lower triangular
318 : !lc> B is a (na,ncb) matrix
319 : !lc> C is a (na,ncb) matrix where optionally only the upper or lower
320 : !lc> triangle may be computed
321 : !lc> \details
322 : !lc> \param uplo_a 'U' if A is upper triangular
323 : !lc> 'L' if A is lower triangular
324 : !lc> anything else if A is a full matrix
325 : !lc> Please note: This pertains to the original A (as set in the calling program)
326 : !lc> whereas the transpose of A is used for calculations
327 : !lc> If uplo_a is 'U' or 'L', the other triangle is not used at all,
328 : !lc> i.e. it may contain arbitrary numbers
329 : !lc> \param uplo_c 'U' if only the upper diagonal part of C is needed
330 : !lc> 'L' if only the upper diagonal part of C is needed
331 : !lc> anything else if the full matrix C is needed
332 : !lc> Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
333 : !lc> written to a certain extent, i.e. one shouldn't rely on the content there!
334 : !lc> \param na Number of rows/columns of A, number of rows of B and C
335 : !lc> \param ncb Number of columns of B and C
336 : !lc> \param a matrix a
337 : !lc> \param lda leading dimension of matrix a
338 : !lc> \param ldaCols columns of matrix a
339 : !lc> \param b matrix b
340 : !lc> \param ldb leading dimension of matrix b
341 : !lc> \param ldbCols columns of matrix b
342 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
343 : !lc> \param mpi_comm_rows MPI communicator for rows
344 : !lc> \param mpi_comm_cols MPI communicator for columns
345 : !lc> \param c matrix c
346 : !lc> \param ldc leading dimension of matrix c
347 : !lc> \result success int report success (1) or failure (0)
348 : !lc> */
349 :
350 : !lc> int elpa_mult_at_b_real_single(char uplo_a, char uplo_c, int na, int ncb, float *a, int lda, int ldaCols, float *b, int ldb, int ldbCols, int nlbk, int mpi_comm_rows, int mpi_comm_cols, float *c, int ldc, int ldcCols);
351 :
352 :
353 : #define REALCASE 1
354 : #define SINGLE_PRECISION 1
355 : #include "../../general/precision_macros.h"
356 : #include "./elpa_mult_at_b_c_interface_template.F90"
357 : #undef SINGLE_PRECISION
358 : #undef REALCASE
359 :
360 : #endif /* WANT_SINGLE_PRECISION_REAL */
361 :
362 : !lc> /*
363 : !lc> \brief C interface for elpa_mult_ah_b_complex_double: Performs C : = A**H * B for double-precision matrices
364 : !lc> where A is a square matrix (na,na) which is optionally upper or lower triangular
365 : !lc> B is a (na,ncb) matrix
366 : !lc> C is a (na,ncb) matrix where optionally only the upper or lower
367 : !lc> triangle may be computed
368 : !lc> \details
369 : !lc>
370 : !lc> \param uplo_a 'U' if A is upper triangular
371 : !lc> 'L' if A is lower triangular
372 : !lc> anything else if A is a full matrix
373 : !lc> Please note: This pertains to the original A (as set in the calling program)
374 : !lc> whereas the transpose of A is used for calculations
375 : !lc> If uplo_a is 'U' or 'L', the other triangle is not used at all,
376 : !lc> i.e. it may contain arbitrary numbers
377 : !lc> \param uplo_c 'U' if only the upper diagonal part of C is needed
378 : !lc> 'L' if only the upper diagonal part of C is needed
379 : !lc> anything else if the full matrix C is needed
380 : !lc> Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
381 : !lc> written to a certain extent, i.e. one shouldn't rely on the content there!
382 : !lc> \param na Number of rows/columns of A, number of rows of B and C
383 : !lc> \param ncb Number of columns of B and C
384 : !lc> \param a matrix a
385 : !lc> \param lda leading dimension of matrix a
386 : !lc> \param b matrix b
387 : !lc> \param ldb leading dimension of matrix b
388 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
389 : !lc> \param mpi_comm_rows MPI communicator for rows
390 : !lc> \param mpi_comm_cols MPI communicator for columns
391 : !lc> \param c matrix c
392 : !lc> \param ldc leading dimension of matrix c
393 : !lc> \result success int reports success (1) or failure (0)
394 : !lc> */
395 :
396 : !lc> int elpa_mult_ah_b_complex_double(char uplo_a, char uplo_c, int na, int ncb, double complex *a, int lda, int ldaCols, double complex *b, int ldb, int ldbCols, int nblk, int mpi_comm_rows, int mpi_comm_cols, double complex *c, int ldc, int ldcCols);
397 : #define COMPLEXCASE 1
398 : #define DOUBLE_PRECISION 1
399 : #include "../../general/precision_macros.h"
400 : #include "./elpa_mult_ah_b_c_interface_template.F90"
401 : #undef DOUBLE_PRECISION
402 : #undef COMPLEXCASE
403 :
404 :
405 : #ifdef WANT_SINGLE_PRECISION_COMPLEX
406 :
407 : !lc> /*
408 : !lc> \brief C interface for elpa_mult_ah_b_complex_single: Performs C : = A**H * B for single-precision matrices
409 : !lc> where A is a square matrix (na,na) which is optionally upper or lower triangular
410 : !lc> B is a (na,ncb) matrix
411 : !lc> C is a (na,ncb) matrix where optionally only the upper or lower
412 : !lc> triangle may be computed
413 : !lc> \details
414 : !lc>
415 : !lc> \param uplo_a 'U' if A is upper triangular
416 : !lc> 'L' if A is lower triangular
417 : !lc> anything else if A is a full matrix
418 : !lc> Please note: This pertains to the original A (as set in the calling program)
419 : !lc> whereas the transpose of A is used for calculations
420 : !lc> If uplo_a is 'U' or 'L', the other triangle is not used at all,
421 : !lc> i.e. it may contain arbitrary numbers
422 : !lc> \param uplo_c 'U' if only the upper diagonal part of C is needed
423 : !lc> 'L' if only the upper diagonal part of C is needed
424 : !lc> anything else if the full matrix C is needed
425 : !lc> Please note: Even when uplo_c is 'U' or 'L', the other triangle may be
426 : !lc> written to a certain extent, i.e. one shouldn't rely on the content there!
427 : !lc> \param na Number of rows/columns of A, number of rows of B and C
428 : !lc> \param ncb Number of columns of B and C
429 : !lc> \param a matrix a
430 : !lc> \param lda leading dimension of matrix a
431 : !lc> \param b matrix b
432 : !lc> \param ldb leading dimension of matrix b
433 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
434 : !lc> \param mpi_comm_rows MPI communicator for rows
435 : !lc> \param mpi_comm_cols MPI communicator for columns
436 : !lc> \param c matrix c
437 : !lc> \param ldc leading dimension of matrix c
438 : !lc> \result success int reports success (1) or failure (0)
439 : !lc> */
440 :
441 : !lc> int elpa_mult_ah_b_complex_single(char uplo_a, char uplo_c, int na, int ncb, complex float *a, int lda, int ldaCols, complex float *b, int ldb, int ldbCols, int nblk, int mpi_comm_rows, int mpi_comm_cols, complex float *c, int ldc, int ldcCols);
442 : #define COMPLEXCASE 1
443 : #define SINGLE_PRECISION 1
444 : #include "../../general/precision_macros.h"
445 : #include "./elpa_mult_ah_b_c_interface_template.F90"
446 : #undef SINGLE_PRECISION
447 : #undef COMPLEXCASE
448 : #endif /* WANT_SINGLE_PRECISION_COMPLEX */
449 :
450 : !lc> /*
451 : !lc> \brief C interface to elpa_invert_trm_real_double: Inverts a real double-precision upper triangular matrix
452 : !lc> \details
453 : !lc> \param na Order of matrix
454 : !lc> \param a(lda,matrixCols) Distributed matrix which should be inverted
455 : !lc> Distribution is like in Scalapack.
456 : !lc> Only upper triangle is needs to be set.
457 : !lc> The lower triangle is not referenced.
458 : !lc> \param lda Leading dimension of a
459 : !lc> \param matrixCols local columns of matrix a
460 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
461 : !lc> \param mpi_comm_rows MPI communicator for rows
462 : !lc> \param mpi_comm_cols MPI communicator for columns
463 : !lc> \param wantDebug int more debug information on failure if 1, else 0
464 : !lc> \result succes int reports success (1) or failure (0)
465 : !lc> */
466 :
467 : !lc> int elpa_invert_trm_real_double(int na, double *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
468 : #define REALCASE 1
469 : #define DOUBLE_PRECISION 1
470 : #include "../../general/precision_macros.h"
471 : #include "./elpa_invert_trm_c_interface_template.F90"
472 : #undef DOUBLE_PRECISION
473 : #undef REALCASE
474 :
475 : #ifdef WANT_SINGLE_PRECISION_REAL
476 :
477 : !lc> /*
478 : !lc> \brief C interface to elpa_invert_trm_real_single: Inverts a real single-precision upper triangular matrix
479 : !lc> \details
480 : !lc> \param na Order of matrix
481 : !lc> \param a(lda,matrixCols) Distributed matrix which should be inverted
482 : !lc> Distribution is like in Scalapack.
483 : !lc> Only upper triangle is needs to be set.
484 : !lc> The lower triangle is not referenced.
485 : !lc> \param lda Leading dimension of a
486 : !lc> \param matrixCols local columns of matrix a
487 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
488 : !lc> \param mpi_comm_rows MPI communicator for rows
489 : !lc> \param mpi_comm_cols MPI communicator for columns
490 : !lc> \param wantDebug int more debug information on failure if 1, else 0
491 : !lc> \result succes int reports success (1) or failure (0)
492 : !lc> */
493 :
494 : !lc> int elpa_invert_trm_real_single(int na, double *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
495 :
496 : #define REALCASE 1
497 : #define SINGLE_PRECISION 1
498 : #include "../../general/precision_macros.h"
499 : #include "./elpa_invert_trm_c_interface_template.F90"
500 : #undef SINGLE_PRECISION
501 : #undef REALCASE
502 :
503 : #endif /* WANT_SINGLE_PRECISION_REAL */
504 :
505 : !lc> /*
506 : !lc> \brief C interface to elpa_invert_trm_complex_double: Inverts a double-precision complex upper triangular matrix
507 : !lc> \details
508 : !lc> \param na Order of matrix
509 : !lc> \param a(lda,matrixCols) Distributed matrix which should be inverted
510 : !lc> Distribution is like in Scalapack.
511 : !lc> Only upper triangle is needs to be set.
512 : !lc> The lower triangle is not referenced.
513 : !lc> \param lda Leading dimension of a
514 : !lc> \param matrixCols local columns of matrix a
515 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
516 : !lc> \param mpi_comm_rows MPI communicator for rows
517 : !lc> \param mpi_comm_cols MPI communicator for columns
518 : !lc> \param wantDebug int more debug information on failure if 1, else 0
519 : !lc> \result succes int reports success (1) or failure (0)
520 : !lc> */
521 :
522 : !lc> int elpa_invert_trm_complex_double(int na, double complex *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
523 : #define COMPLEXCASE 1
524 : #define DOUBLE_PRECISION 1
525 : #include "../../general/precision_macros.h"
526 : #include "./elpa_invert_trm_c_interface_template.F90"
527 : #undef DOUBLE_PRECISION
528 : #undef COMPLEXCASE
529 :
530 : #ifdef WANT_SINGLE_PRECISION_COMPLEX
531 : !lc> /*
532 : !lc> \brief C interface to elpa_invert_trm_complex_single: Inverts a single-precision complex upper triangular matrix
533 : !lc> \details
534 : !lc> \param na Order of matrix
535 : !lc> \param a(lda,matrixCols) Distributed matrix which should be inverted
536 : !lc> Distribution is like in Scalapack.
537 : !lc> Only upper triangle is needs to be set.
538 : !lc> The lower triangle is not referenced.
539 : !lc> \param lda Leading dimension of a
540 : !lc> \param matrixCols local columns of matrix a
541 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
542 : !lc> \param mpi_comm_rows MPI communicator for rows
543 : !lc> \param mpi_comm_cols MPI communicator for columns
544 : !lc> \param wantDebug int more debug information on failure if 1, else 0
545 : !lc> \result succes int reports success (1) or failure (0)
546 : !lc> */
547 :
548 : !lc> int elpa_invert_trm_complex_single(int na, complex float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
549 : #define COMPLEXCASE 1
550 : #define SINGLE_PRECISION 1
551 : #include "../../general/precision_macros.h"
552 : #include "./elpa_invert_trm_c_interface_template.F90"
553 : #undef SINGLE_PRECISION
554 : #undef COMPLEXCASE
555 :
556 :
557 : #endif /* WANT_SINGLE_PRECISION_COMPLEX */
558 :
559 : !lc> /*
560 : !lc> \brief elpa_cholesky_real_double: Cholesky factorization of a double-precision real symmetric matrix
561 : !lc> \details
562 : !lc>
563 : !lc> *\param na Order of matrix
564 : !lc> *\param a(lda,matrixCols) Distributed matrix which should be factorized.
565 : !lc> * Distribution is like in Scalapack.
566 : !lc> * Only upper triangle is needs to be set.
567 : !lc> * On return, the upper triangle contains the Cholesky factor
568 : !lc> * and the lower triangle is set to 0.
569 : !lc> *\param lda Leading dimension of a
570 : !lc> *\param matrixCols local columns of matrix a
571 : !lc> *\param nblk blocksize of cyclic distribution, must be the same in both directions!
572 : !lc> *\param mpi_comm_rows MPI communicator for rows
573 : !lc> *\param mpi_comm_cols MPI communicator for columns
574 : !lc> *\param wantDebug int more debug information on failure if 1, else 0
575 : !lc> *\result succes int reports success (1) or failure (0)
576 : !lc> */
577 :
578 : !lc> int elpa_cholesky_real_double(int na, double *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
579 : #define REALCASE 1
580 : #define DOUBLE_PRECISION 1
581 : #include "../../general/precision_macros.h"
582 : #include "./elpa_cholesky_c_interface_template.F90"
583 : #undef DOUBLE_PRECISION
584 : #undef REALCASE
585 :
586 : #ifdef WANT_SINGLE_PRECISION_REAL
587 :
588 : !lc> /*
589 : !lc> \brief elpa_cholesky_real_single: Cholesky factorization of a single-precision real symmetric matrix
590 : !lc> \details
591 : !lc>
592 : !lc> \param na Order of matrix
593 : !lc> \param a(lda,matrixCols) Distributed matrix which should be factorized.
594 : !lc> Distribution is like in Scalapack.
595 : !lc> Only upper triangle is needs to be set.
596 : !lc> On return, the upper triangle contains the Cholesky factor
597 : !lc> and the lower triangle is set to 0.
598 : !lc> \param lda Leading dimension of a
599 : !lc> \param matrixCols local columns of matrix a
600 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
601 : !lc> \param mpi_comm_rows MPI communicator for rows
602 : !lc> \param mpi_comm_cols MPI communicator for columns
603 : !lc> \param wantDebug int more debug information on failure if 1, else 0
604 : !lc> \result succes int reports success (1) or failure (0)
605 : !lc> */
606 :
607 : !lc> int elpa_cholesky_real_single(int na, float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
608 : #define REALCASE 1
609 : #define SINGLE_PRECISION 1
610 : #include "../../general/precision_macros.h"
611 : #include "./elpa_cholesky_c_interface_template.F90"
612 : #undef SINGLE_PRECISION
613 : #undef REALCASE
614 :
615 :
616 : #endif /* WANT_SINGLE_PRECISION_REAL */
617 :
618 : !lc> /*
619 : !lc> \brief C interface elpa_cholesky_complex_double: Cholesky factorization of a double-precision complex hermitian matrix
620 : !lc> \details
621 : !lc> \param na Order of matrix
622 : !lc> \param a(lda,matrixCols) Distributed matrix which should be factorized.
623 : !lc> Distribution is like in Scalapack.
624 : !lc> Only upper triangle is needs to be set.
625 : !lc> On return, the upper triangle contains the Cholesky factor
626 : !lc> and the lower triangle is set to 0.
627 : !lc> \param lda Leading dimension of a
628 : !lc> \param matrixCols local columns of matrix a
629 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
630 : !lc> \param mpi_comm_rows MPI communicator for rows
631 : !lc> \param mpi_comm_cols MPI communicator for columns
632 : !lc> \param wantDebug int more debug information on failure, if 1, else 0
633 : !lc> \result succes int reports success (1) or failure (0)
634 : !lc> */
635 :
636 : !lc> int elpa_cholesky_complex_double(int na, double complex *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
637 : #define COMPLEXCASE 1
638 : #define DOUBLE_PRECISION 1
639 : #include "../../general/precision_macros.h"
640 : #include "./elpa_cholesky_c_interface_template.F90"
641 : #undef DOUBLE_PRECISION
642 : #undef COMPLEXCASE
643 :
644 : #ifdef WANT_SINGLE_PRECISION_COMPLEX
645 :
646 : !lc> /*
647 : !lc> \brief C interface elpa_cholesky_complex_single: Cholesky factorization of a single-precision complex hermitian matrix
648 : !lc> \details
649 : !lc> \param na Order of matrix
650 : !lc> \param a(lda,matrixCols) Distributed matrix which should be factorized.
651 : !lc> Distribution is like in Scalapack.
652 : !lc> Only upper triangle is needs to be set.
653 : !lc> On return, the upper triangle contains the Cholesky factor
654 : !lc> and the lower triangle is set to 0.
655 : !lc> \param lda Leading dimension of a
656 : !lc> \param matrixCols local columns of matrix a
657 : !lc> \param nblk blocksize of cyclic distribution, must be the same in both directions!
658 : !lc> \param mpi_comm_rows MPI communicator for rows
659 : !lc> \param mpi_comm_cols MPI communicator for columns
660 : !lc> \param wantDebug int more debug information on failure, if 1, else 0
661 : !lc> \result succes int reports success (1) or failure (0)
662 : !lc> */
663 :
664 : !lc> int elpa_cholesky_complex_single(int na, complex float *a, int lda, int nblk, int matrixCols, int mpi_comm_rows, int mpi_comm_cols, int wantDebug);
665 : #define COMPLEXCASE 1
666 : #define SINGLE_PRECISION 1
667 : #include "../../general/precision_macros.h"
668 : #include "./elpa_cholesky_c_interface_template.F90"
669 : #undef SINGLE_PRECISION
670 : #undef COMPLEXCASE
671 : #endif /* WANT_SINGLE_PRECISION_COMPLEX */
|