KIDS  ver-0.0.1
KIDS : Kernel Integrated Dynamics Simulator
Loading...
Searching...
No Matches
linalg.h
Go to the documentation of this file.
1
30#ifndef KIDS_LINALG_H
31#define KIDS_LINALG_H
32
33#include "kids/Types.h"
34
35#define KIDS_LINALG_BIND_EIGEN_NOT_USE_TEMPLATE
36
37#ifndef KIDS_LINALG_BIND_EIGEN_NOT_USE_TEMPLATE
38#include "kids/linalg_tpl.h"
39#else // KIDS_LINALG_BIND_EIGEN_NOT_USE_TEMPLATE
40
41namespace PROJECT_NS {
42
50bool ARRAY_ISFINITE(kids_real* A, size_t n);
51
59bool ARRAY_ISFINITE(kids_complex* A, size_t n);
60
67void ARRAY_CLEAR(kids_int* A, size_t N);
68
75void ARRAY_CLEAR(kids_real* A, size_t N);
76
83void ARRAY_CLEAR(kids_complex* A, size_t N);
84
95void ARRAY_MATMUL(kids_real* A, kids_real* B, kids_real* C, size_t N1, size_t N2, size_t N3);
96
107void ARRAY_MATMUL(kids_complex* A, kids_complex* B, kids_complex* C, size_t N1, size_t N2, size_t N3);
108
119void ARRAY_MATMUL(kids_complex* A, kids_real* B, kids_complex* C, size_t N1, size_t N2, size_t N3);
120
131void ARRAY_MATMUL(kids_complex* A, kids_complex* B, kids_real* C, size_t N1, size_t N2, size_t N3);
132
143void ARRAY_MATMUL_TRANS1(kids_real* A, kids_real* B, kids_real* C, size_t N1, size_t N2, size_t N3);
144
155void ARRAY_MATMUL_TRANS1(kids_complex* A, kids_complex* B, kids_complex* C, size_t N1, size_t N2, size_t N3);
156
167void ARRAY_MATMUL_TRANS1(kids_complex* A, kids_real* B, kids_complex* C, size_t N1, size_t N2, size_t N3);
168
179void ARRAY_MATMUL_TRANS1(kids_complex* A, kids_complex* B, kids_real* C, size_t N1, size_t N2, size_t N3);
180
191void ARRAY_MATMUL_TRANS2(kids_real* A, kids_real* B, kids_real* C, size_t N1, size_t N2, size_t N3);
192
203void ARRAY_MATMUL_TRANS2(kids_complex* A, kids_complex* B, kids_complex* C, size_t N1, size_t N2, size_t N3);
204
215void ARRAY_MATMUL_TRANS2(kids_complex* A, kids_real* B, kids_complex* C, size_t N1, size_t N2, size_t N3);
216
227void ARRAY_MATMUL_TRANS2(kids_complex* A, kids_complex* B, kids_real* C, size_t N1, size_t N2, size_t N3);
228
238void ARRAY_OUTER_TRANS2(kids_real* A, kids_real* B, kids_real* C, size_t N1, size_t N2);
239
249void ARRAY_OUTER_TRANS2(kids_complex* A, kids_complex* B, kids_complex* C, size_t N1, size_t N2);
250
264 size_t N1, size_t N2, size_t N0, size_t N3);
265
279 size_t N1, size_t N2, size_t N0, size_t N3);
280
294 size_t N1, size_t N2, size_t N0, size_t N3);
295
309 size_t N1, size_t N2, size_t N0, size_t N3);
310
324 size_t N1, size_t N2, size_t N0, size_t N3);
325
339 size_t N1, size_t N2, size_t N0, size_t N3);
340
354 size_t N1, size_t N2, size_t N0, size_t N3);
355
369 size_t N1, size_t N2, size_t N0, size_t N3);
370
380kids_real ARRAY_TRACE2(kids_real* B, kids_real* C, size_t N1, size_t N2);
381
391kids_complex ARRAY_TRACE2(kids_complex* B, kids_complex* C, size_t N1, size_t N2);
392
402kids_complex ARRAY_TRACE2(kids_complex* B, kids_real* C, size_t N1, size_t N2);
403
413kids_complex ARRAY_TRACE2(kids_real* B, kids_complex* C, size_t N1, size_t N2);
414
424kids_real ARRAY_TRACE2_DIAG(kids_real* B, kids_real* C, size_t N1, size_t N2);
425
435kids_complex ARRAY_TRACE2_DIAG(kids_complex* B, kids_complex* C, size_t N1, size_t N2);
436
446kids_complex ARRAY_TRACE2_DIAG(kids_complex* B, kids_real* C, size_t N1, size_t N2);
447
457kids_complex ARRAY_TRACE2_DIAG(kids_real* B, kids_complex* C, size_t N1, size_t N2);
458
468kids_real ARRAY_TRACE2_OFFD(kids_real* B, kids_real* C, size_t N1, size_t N2);
469
479kids_complex ARRAY_TRACE2_OFFD(kids_complex* B, kids_complex* C, size_t N1, size_t N2);
480
490kids_complex ARRAY_TRACE2_OFFD(kids_complex* B, kids_real* C, size_t N1, size_t N2);
491
501kids_complex ARRAY_TRACE2_OFFD(kids_real* B, kids_complex* C, size_t N1, size_t N2);
502
512
522
532
542
554kids_real ARRAY_INNER_VMV_TRANS1(kids_real* B, kids_real* C, kids_real* D, size_t N1, size_t N2);
555
568
581
593kids_complex ARRAY_INNER_VMV_TRANS1(kids_real* B, kids_complex* C, kids_real* D, size_t N1, size_t N2);
594
601void ARRAY_EYE(kids_real* A, size_t n);
602
609void ARRAY_EYE(kids_complex* A, size_t n);
610
618void ARRAY_MAT_DIAG(kids_real* A, kids_real* B, size_t N1);
619
627void ARRAY_MAT_DIAG(kids_complex* A, kids_complex* B, size_t N1);
628
636void ARRAY_MAT_DIAG(kids_complex* A, kids_real* B, size_t N1);
637
645void ARRAY_MAT_OFFD(kids_real* A, kids_real* B, size_t N1);
646
654void ARRAY_MAT_OFFD(kids_complex* A, kids_complex* B, size_t N1);
655
663void ARRAY_MAT_OFFD(kids_complex* A, kids_real* B, size_t N1);
664
673void LinearSolve(kids_real* x, kids_real* A, kids_real* b, size_t N);
674
683void EigenSolve(kids_real* E, kids_real* T, kids_real* A, size_t N);
684
693void EigenSolve(kids_real* E, kids_complex* T, kids_complex* A, size_t N);
694
703void EigenSolve(kids_complex* E, kids_complex* T, kids_complex* A, size_t N);
704
713void PseudoInverse(kids_real* A, kids_real* invA, size_t N, kids_real e = 1E-5);
714
722void MatrixInverse(kids_real* invA, kids_real* A, size_t N);
723
731void ARRAY_INV_MAT(kids_complex* invA, kids_complex* A, size_t N);
732
733void ARRAY_EXP_MAT_GENERAL(kids_complex* expkA, kids_complex* A, kids_complex k, size_t N);
734
735void ARRAY_CORRECT_U(kids_complex* U, size_t N);
736
737void ARRAY_TRANSPOSE(kids_real* A, size_t N1, size_t N2);
738
739void ARRAY_TRANSPOSE(kids_complex* A, size_t N1, size_t N2);
740
741
742}; // namespace PROJECT_NS
743
744#endif // KIDS_LINALG_BIND_EIGEN_NOT_USE_TEMPLATE
745
746#endif // KIDS_LINALG_H
definition of types in the project and some utiles for types
#define ARRAY_MATMUL_TRANS2(_A, _B, _C, _n1, _n2, _n3)
#define ARRAY_MATMUL3_TRANS1(_A, _B, _C, _D, _n1, _n2, _n0, _n3)
#define ARRAY_MATMUL3_TRANS2(_A, _B, _C, _D, _n1, _n2, _n0, _n3)
#define ARRAY_MATMUL(_A, _B, _C, _n1, _n2, _n3)
#define ARRAY_EYE(_A, _n)
Definition array_macro.h:75
#define ARRAY_MATMUL_TRANS1(_A, _B, _C, _n1, _n2, _n3)
#define ARRAY_TRACE2(_B, _C, _n1, _n2)
#define ARRAY_CLEAR(_A, _n)
Definition array_macro.h:36
Provide linalg APIs by templating.
< http://warp.povusers.org/FunctionParser/fparser.html
Definition Context.h:39
bool ARRAY_ISFINITE(kids_real *A, size_t n)
Check if all elements of a real array are finite.
Definition linalg.cpp:38
void PseudoInverse(kids_real *A, kids_real *invA, size_t N, kids_real e=1E-5)
Compute the pseudo-inverse of a matrix for real numbers.
Definition linalg_1.cpp:71
void ARRAY_CORRECT_U(kids_complex *U, size_t N)
Definition linalg_1.cpp:132
void ARRAY_MAT_DIAG(kids_real *A, kids_real *B, size_t N1)
Copy the diagonal elements from matrix B to matrix A for real matrices.
Definition linalg.cpp:458
kids_real ARRAY_INNER_TRANS1(kids_real *B, kids_real *C, size_t N1)
Compute the inner product of the transpose of B with C for real matrices.
Definition linalg.cpp:400
void ARRAY_INV_MAT(kids_complex *invA, kids_complex *A, size_t N)
Compute the inverse of a matrix for complex numbers.
Definition linalg_1.cpp:91
void ARRAY_TRANSPOSE(kids_real *A, size_t N1, size_t N2)
Definition linalg.cpp:500
void ARRAY_MAT_OFFD(kids_real *A, kids_real *B, size_t N1)
Copy the off-diagonal elements from matrix B to matrix A for real matrices.
Definition linalg.cpp:479
double kids_real
Alias for real number type.
Definition Types.h:59
std::complex< double > kids_complex
Alias for complex number type.
Definition Types.h:60
void ARRAY_OUTER_TRANS2(kids_real *A, kids_real *B, kids_real *C, size_t N1, size_t N2)
Perform outer product with transpose where A = B^T * C for real matrices.
Definition linalg.cpp:160
kids_real ARRAY_INNER_VMV_TRANS1(kids_real *B, kids_real *C, kids_real *D, size_t N1, size_t N2)
Compute the inner product of the transpose of the real vector B with the matrix C,...
Definition linalg.cpp:424
void EigenSolve(kids_real *E, kids_real *T, kids_real *A, size_t N)
Solve the eigenvalue problem for real matrices.
Definition linalg_1.cpp:44
kids_real ARRAY_TRACE2_DIAG(kids_real *B, kids_real *C, size_t N1, size_t N2)
Compute the trace of the diagonal elements of the matrix product B * C for real matrices.
Definition linalg.cpp:340
void ARRAY_EXP_MAT_GENERAL(kids_complex *expkA, kids_complex *A, kids_complex k, size_t N)
Definition linalg_1.cpp:97
void MatrixInverse(kids_real *invA, kids_real *A, size_t N)
Compute the inverse of a matrix for real numbers.
void LinearSolve(kids_real *x, kids_real *A, kids_real *b, size_t N)
Solve a linear system Ax = b for real matrices.
Definition linalg_1.cpp:37
kids_real ARRAY_TRACE2_OFFD(kids_real *B, kids_real *C, size_t N1, size_t N2)
Compute the trace of the off-diagonal elements of the matrix product B * C for real matrices.
Definition linalg.cpp:368
int kids_int
Alias for integer type.
Definition Types.h:57