30#ifndef KIDS_LINALG_TPL_H
31#define KIDS_LINALG_TPL_H
42#define EigMajor Eigen::RowMajor
45using EigVX = Eigen::Matrix<T, Eigen::Dynamic, 1, EigMajor>;
48using EigMX = Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, EigMajor>;
51using EigAX = Eigen::Array<T, Eigen::Dynamic, Eigen::Dynamic, EigMajor>;
68 for (
int i = 0; i < n; ++i)
69 if (!std::isfinite(std::abs(A[i])))
return false;
75 memset(A, 0, N *
sizeof(T));
78template <
class TA,
class TB,
class TC>
79void ARRAY_MATMUL(TA* A, TB* B, TC* C,
size_t N1,
size_t N2,
size_t N3) {
80 Eigen::Map<EigMX<TA>> MapA(A, N1, N3);
81 Eigen::Map<EigMX<TB>> MapB(B, N1, N2);
82 Eigen::Map<EigMX<TC>> MapC(C, N2, N3);
83 MapA = (MapB * MapC).eval();
86template <
class TA,
class TB,
class TC>
88 Eigen::Map<EigMX<TA>> MapA(A, N1, N3);
89 Eigen::Map<EigMX<TB>> MapB(B, N2, N1);
90 Eigen::Map<EigMX<TC>> MapC(C, N2, N3);
91 MapA = (MapB.adjoint() * MapC).eval();
94template <
class TA,
class TB,
class TC>
96 Eigen::Map<EigMX<TA>> MapA(A, N1, N3);
97 Eigen::Map<EigMX<TB>> MapB(B, N1, N2);
98 Eigen::Map<EigMX<TC>> MapC(C, N3, N2);
99 MapA = (MapB * MapC.adjoint()).eval();
107template <
class TA,
class T,
class TC>
110 Eigen::Map<EigMX<TA>> MapA(A, N1, N3);
111 Eigen::Map<EigMX<T>> MapB(B, N2, N1);
112 Eigen::Map<EigMX<TC>> MapC(C, N2, 1);
113 Eigen::Map<EigMX<T>> MapD(D, N2, N3);
114 MapA = (MapB.adjoint() * (MapC.asDiagonal() * MapD)).eval();
116 Eigen::Map<EigMX<TA>> MapA(A, N1, N3);
117 Eigen::Map<EigMX<T>> MapB(B, N2, N1);
118 Eigen::Map<EigMX<TC>> MapC(C, N2, N2);
119 Eigen::Map<EigMX<T>> MapD(D, N2, N3);
120 MapA = (MapB.adjoint() * (MapC * MapD)).eval();
124template <
class TA,
class T,
class TC>
127 Eigen::Map<EigMX<TA>> MapA(A, N1, N3);
128 Eigen::Map<EigMX<T>> MapB(B, N1, N2);
129 Eigen::Map<EigMX<TC>> MapC(C, N2, 1);
130 Eigen::Map<EigMX<T>> MapD(D, N3, N2);
131 MapA = (MapB * (MapC.asDiagonal() * MapD.adjoint())).eval();
133 Eigen::Map<EigMX<TA>> MapA(A, N1, N3);
134 Eigen::Map<EigMX<T>> MapB(B, N1, N2);
135 Eigen::Map<EigMX<TC>> MapC(C, N2, N2);
136 Eigen::Map<EigMX<T>> MapD(D, N3, N2);
137 MapA = (MapB * (MapC * MapD.adjoint())).eval();
141template <
class TB,
class TC>
143 Eigen::Map<EigMX<TB>> MapB(B, N1, N2);
144 Eigen::Map<EigMX<TC>> MapC(C, N2, N1);
145 TB res = (MapB.array() * (MapC.transpose()).array()).sum();
149template <
class TB,
class TC>
151 Eigen::Map<EigMX<TB>> MapB(B, N1, N2);
152 Eigen::Map<EigMX<TC>> MapC(C, N2, N1);
153 TB res = (MapB.diagonal().array() * MapC.diagonal().array()).sum();
157template <
class TB,
class TC>
159 Eigen::Map<EigMX<TB>> MapB(B, N1, N2);
160 Eigen::Map<EigMX<TC>> MapC(C, N2, N1);
161 TB res = (MapB.array() * (MapC.transpose()).array()).sum()
162 - (MapB.diagonal().array() * MapC.diagonal().array()).sum();
166template <
class TB,
class TC>
168 Eigen::Map<EigMX<TB>> MapB(B, N1, 1);
169 Eigen::Map<EigMX<TC>> MapC(C, N1, 1);
170 return (MapB.adjoint() * MapC).sum();
173template <
class TB,
class TC,
class TD>
175 Eigen::Map<EigMX<TB>> MapB(B, N1, 1);
176 Eigen::Map<EigMX<TC>> MapC(C, N1, N2);
177 Eigen::Map<EigMX<TD>> MapD(D, N2, 1);
178 return (MapB.adjoint() * MapC * MapD).sum();
183 Eigen::Map<EigMX<T>> MapA(A, n, n);
187template <
class TA,
class TB>
189 Eigen::Map<EigMX<TA>> MapA(A, N1, N1);
190 Eigen::Map<EigMX<TB>> MapB(B, N1, N1);
193 MapA.diagonal() = MapB.diagonal();
196template <
class TA,
class TB>
198 Eigen::Map<EigMX<TA>> MapA(A, N1, N1);
199 Eigen::Map<EigMX<TB>> MapB(B, N1, N1);
201 MapA.diagonal().array() = TA(0);
212 Mapx = MapA.householderQr().solve(Mapb);
219 Eigen::SelfAdjointEigenSolver<EigMXr> eig(MapA);
220 MapE = eig.eigenvalues().real();
221 MapT = eig.eigenvectors().real();
228 Eigen::SelfAdjointEigenSolver<EigMXc> eig(MapA);
229 MapE = eig.eigenvalues().real();
230 MapT = eig.eigenvectors();
237 Eigen::ComplexEigenSolver<EigMXc> eig(MapA);
238 MapE = eig.eigenvalues();
239 MapT = eig.eigenvectors();
244 MapMXr MapInvA(invA, N, N);
245 MapInvA = MapA.completeOrthogonalDecomposition().pseudoInverse();
420A.cast<kids_real>(); // kids_real(A)
421A.cast<float>(); // single(A)
422A.cast<int>(); // int32(A)
425// if the original type equals destination type, no work is done
428// Solve Ax = b. Result stored in x. Matlab: x = A \ b.
429x = A.ldlt().solve(b)); // A sym. p.s.d. #include <Eigen/Cholesky>
430x = A.llt() .solve(b)); // A sym. p.d. #include <Eigen/Cholesky>
431x = A.lu() .solve(b)); // Stable and fast. #include <Eigen/LU>
432x = A.qr() .solve(b)); // No pivoting. #include <Eigen/QR>
433x = A.svd() .solve(b)); // Stable, slowest. #include <Eigen/SVD>
434// .ldlt() -> .matrixL() and .matrixD()
435// .llt() -> .matrixL()
436// .lu() -> .matrixL() and .matrixU()
437// .qr() -> .matrixQ() and .matrixR()
438// .svd() -> .matrixU(), .singularValues(), and .matrixV()
440// Eigenvalue problems
442A.eigenvalues(); // eig(A);
443EigenSolver<Matrix3d> eig(A); // [vec val] = eig(A)
444eig.eigenvalues(); // diag(val)
445eig.eigenvectors(); // vec
446// For self-adjoint matrices use SelfAdjointEigenSolver<>
449// Use Map of Eigen, associate Map class with C++ pointer
452std::cout << std::setiosflags(std::ios::scientific)
453 << std::setprecision(8) << W << std::endl;
455*****************************************/
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)
#define ARRAY_MATMUL_TRANS1(_A, _B, _C, _n1, _n2, _n3)
#define ARRAY_TRACE2(_B, _C, _n1, _n2)
#define ARRAY_CLEAR(_A, _n)
EigMX< kids_real > EigMXr
EigMX< kids_complex > EigAXc
Eigen::Map< EigMXr > MapMXr
Eigen::Map< EigMXc > MapMXc
Eigen::Map< EigAXr > MapAXr
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, EigMajor > EigMX
Eigen::Map< EigVXr > MapVXr
EigMX< kids_real > EigAXr
Eigen::Map< EigAXc > MapAXc
EigVX< kids_complex > EigVXc
Eigen::Matrix< T, Eigen::Dynamic, 1, EigMajor > EigVX
Eigen::Array< T, Eigen::Dynamic, Eigen::Dynamic, EigMajor > EigAX
EigVX< kids_real > EigVXr
EigMX< kids_complex > EigMXc
Eigen::Map< EigVXc > MapVXc
< http://warp.povusers.org/FunctionParser/fparser.html
bool ARRAY_ISFINITE(kids_real *A, size_t n)
Check if all elements of a real array are finite.
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.
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.
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.
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.
double kids_real
Alias for real number type.
std::complex< double > kids_complex
Alias for complex number type.
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.
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,...
void EigenSolve(kids_real *E, kids_real *T, kids_real *A, size_t N)
Solve the eigenvalue problem for real matrices.
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.
void LinearSolve(kids_real *x, kids_real *A, kids_real *b, size_t N)
Solve a linear system Ax = b for real matrices.
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.