KIDS  ver-0.0.1
KIDS : Kernel Integrated Dynamics Simulator
Loading...
Searching...
No Matches
linalg_tpl.h
Go to the documentation of this file.
1
30#ifndef KIDS_LINALG_TPL_H
31#define KIDS_LINALG_TPL_H
32
33#include <cmath>
34#include <complex>
35
36#include "Eigen/Dense"
37#include "Eigen/QR"
38#include "kids/Types.h"
39
40using namespace PROJECT_NS;
41
42#define EigMajor Eigen::RowMajor
43
44template <class T>
45using EigVX = Eigen::Matrix<T, Eigen::Dynamic, 1, EigMajor>;
46
47template <class T>
48using EigMX = Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, EigMajor>;
49
50template <class T>
51using EigAX = Eigen::Array<T, Eigen::Dynamic, Eigen::Dynamic, EigMajor>;
52
59typedef Eigen::Map<EigVXr> MapVXr;
60typedef Eigen::Map<EigVXc> MapVXc;
61typedef Eigen::Map<EigMXr> MapMXr;
62typedef Eigen::Map<EigMXc> MapMXc;
63typedef Eigen::Map<EigAXr> MapAXr;
64typedef Eigen::Map<EigAXc> MapAXc;
65
66template <class T>
67bool ARRAY_ISFINITE(T* A, size_t n) {
68 for (int i = 0; i < n; ++i)
69 if (!std::isfinite(std::abs(A[i]))) return false;
70 return true;
71}
72
73template <class T>
74void ARRAY_CLEAR(T* A, size_t N) {
75 memset(A, 0, N * sizeof(T));
76}
77
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();
84}
85
86template <class TA, class TB, class TC>
87void ARRAY_MATMUL_TRANS1(TA* A, TB* B, TC* C, size_t N1, size_t N2, size_t N3) {
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();
92}
93
94template <class TA, class TB, class TC>
95void ARRAY_MATMUL_TRANS2(TA* A, TB* B, TC* C, size_t N1, size_t N2, size_t N3) {
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();
100}
101
102template <class T>
103void ARRAY_OUTER_TRANS2(T* A, T* B, T* C, size_t N1, size_t N2) {
104 ARRAY_MATMUL_TRANS2(A, B, C, N1, 1, N2);
105}
106
107template <class TA, class T, class TC>
108void ARRAY_MATMUL3_TRANS1(TA* A, T* B, TC* C, T* D, size_t N1, size_t N2, size_t N0, size_t N3) {
109 if (N0 == 0) {
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();
115 } else { // N0 == N2
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();
121 }
122}
123
124template <class TA, class T, class TC>
125void ARRAY_MATMUL3_TRANS2(TA* A, T* B, TC* C, T* D, size_t N1, size_t N2, size_t N0, size_t N3) {
126 if (N0 == 0) {
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();
132 } else { // N0 == N2
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();
138 }
139}
140
141template <class TB, class TC>
142TB ARRAY_TRACE2(TB* B, TC* C, size_t N1, size_t N2) {
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();
146 return res;
147}
148
149template <class TB, class TC>
150TB ARRAY_TRACE2_DIAG(TB* B, TC* C, size_t N1, size_t N2) {
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();
154 return res;
155}
156
157template <class TB, class TC>
158TB ARRAY_TRACE2_OFFD(TB* B, TC* C, size_t N1, size_t N2) {
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();
163 return res;
164}
165
166template <class TB, class TC>
167TB ARRAY_INNER_TRANS1(TB* B, TC* C, size_t N1) {
168 Eigen::Map<EigMX<TB>> MapB(B, N1, 1);
169 Eigen::Map<EigMX<TC>> MapC(C, N1, 1);
170 return (MapB.adjoint() * MapC).sum();
171}
172
173template <class TB, class TC, class TD>
174TB ARRAY_INNER_VMV_TRANS1(TB* B, TC* C, TD* D, size_t N1, size_t N2) {
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();
179}
180
181template <class T>
182void ARRAY_EYE(T* A, size_t n) {
183 Eigen::Map<EigMX<T>> MapA(A, n, n);
184 MapA = EigMX<T>::Identity(n, n);
185}
186
187template <class TA, class TB>
188void ARRAY_MAT_DIAG(TA* A, TB* B, size_t N1) {
189 Eigen::Map<EigMX<TA>> MapA(A, N1, N1);
190 Eigen::Map<EigMX<TB>> MapB(B, N1, N1);
191 ARRAY_CLEAR(A, N1 * N1);
192 // MapA.array() = TA(0);
193 MapA.diagonal() = MapB.diagonal();
194}
195
196template <class TA, class TB>
197void ARRAY_MAT_OFFD(TA* A, TB* B, size_t N1) {
198 Eigen::Map<EigMX<TA>> MapA(A, N1, N1);
199 Eigen::Map<EigMX<TB>> MapB(B, N1, N1);
200 MapA = MapB;
201 MapA.diagonal().array() = TA(0);
202}
203
204/*===========================================================
205= realize interface of linear algebra =
206===========================================================*/
207
208inline void LinearSolve(kids_real* x, kids_real* A, kids_real* b, size_t N) {
209 MapMXr Mapx(x, N, 1);
210 MapMXr MapA(A, N, N);
211 MapMXr Mapb(b, N, 1);
212 Mapx = MapA.householderQr().solve(Mapb);
213}
214
215inline void EigenSolve(kids_real* E, kids_real* T, kids_real* A, size_t N) {
216 MapMXr MapE(E, N, 1);
217 MapMXr MapT(T, N, N);
218 MapMXr MapA(A, N, N);
219 Eigen::SelfAdjointEigenSolver<EigMXr> eig(MapA);
220 MapE = eig.eigenvalues().real();
221 MapT = eig.eigenvectors().real();
222}
223
224inline void EigenSolve(kids_real* E, kids_complex* T, kids_complex* A, size_t N) {
225 MapMXr MapE(E, N, 1);
226 MapMXc MapT(T, N, N);
227 MapMXc MapA(A, N, N);
228 Eigen::SelfAdjointEigenSolver<EigMXc> eig(MapA);
229 MapE = eig.eigenvalues().real();
230 MapT = eig.eigenvectors();
231}
232
233inline void EigenSolve(kids_complex* E, kids_complex* T, kids_complex* A, size_t N) {
234 MapMXc MapE(E, N, 1);
235 MapMXc MapT(T, N, N);
236 MapMXc MapA(A, N, N);
237 Eigen::ComplexEigenSolver<EigMXc> eig(MapA);
238 MapE = eig.eigenvalues();
239 MapT = eig.eigenvectors();
240}
241
242inline void PseudoInverse(kids_real* A, kids_real* invA, size_t N, kids_real e = 1E-5) {
243 MapMXr MapA(A, N, N);
244 MapMXr MapInvA(invA, N, N);
245 MapInvA = MapA.completeOrthogonalDecomposition().pseudoInverse();
246 /* Eigen::JacobiSVD<EigMXr> svd = MapA.jacobiSvd(Eigen::ComputeFullU|Eigen::ComputeFullV);
247 EigMXr invS = EigMXr::Zero(N, N);
248 for (int i = 0; i < N; ++i) {
249 if (svd.singularValues()(i) > e || svd.singularValues()(i) < -e) {
250 invS(i, i) = 1 / svd.singularValues()(i);
251 }
252 }
253 MapInvA = svd.matrixV()*invS*svd.matrixU().transpose(); */
254}
255
256/*****************************************
257 An brief introduction to Eigen library
258------------------------------------------
259
260// Eigen // Matlab // comments
261x.size() // length(x) // vector size
262C.rows() // size(C,1) // number of rows
263C.cols() // size(C,2) // number of columns
264x(i) // x(i+1) // Matlab is 1-based
265C(i, j) // C(i+1,j+1) //
266
267A.resize(4, 4); // Runtime error if assertions are on.
268B.resize(4, 9); // Runtime error if assertions are on.
269A.resize(3, 3); // Ok; size didn't change.
270B.resize(3, 9); // Ok; only dynamic cols changed.
271
272A << 1, 2, 3, // Initialize A. The elements can also be
273 4, 5, 6, // matrices, which are stacked along cols
274 7, 8, 9; // and then the rows are stacked.
275B << A, A, A; // B is three horizontally stacked A's.
276A.fill(10); // Fill A with all 10's.
277
278
279// Eigen // Matlab
280MatrixXd::Identity(rows,cols) // eye(rows,cols)
281C.setIdentity(rows,cols) // C = eye(rows,cols)
282MatrixXd::Zero(rows,cols) // zeros(rows,cols)
283C.setZero(rows,cols) // C = ones(rows,cols)
284MatrixXd::Ones(rows,cols) // ones(rows,cols)
285C.setOnes(rows,cols) // C = ones(rows,cols)
286MatrixXd::Random(rows,cols) // rand(rows,cols)*2-1
287C.setRandom(rows,cols) // C = rand(rows,cols)*2-1
288VectorXd::LinSpaced(size,low,high) // linspace(low,high,size)'
289v.setLinSpaced(size,low,high) // v = linspace(low,high,size)'
290
291
292// Matrix slicing and blocks. All expressions listed here are read/write.
293// Templated size versions are faster. Note that Matlab is 1-based (a size N
294// vector is x(1)...x(N)).
295// Eigen // Matlab
296x.head(n) // x(1:n)
297x.head<n>() // x(1:n)
298x.tail(n) // x(end - n + 1: end)
299x.tail<n>() // x(end - n + 1: end)
300x.segment(i, n) // x(i+1 : i+n)
301x.segment<n>(i) // x(i+1 : i+n)
302P.block(i, j, rows, cols) // P(i+1 : i+rows, j+1 : j+cols)
303P.block<rows, cols>(i, j) // P(i+1 : i+rows, j+1 : j+cols)
304P.row(i) // P(i+1, :)
305P.col(j) // P(:, j+1)
306P.leftCols<cols>() // P(:, 1:cols)
307P.leftCols(cols) // P(:, 1:cols)
308P.middleCols<cols>(j) // P(:, j+1:j+cols)
309P.middleCols(j, cols) // P(:, j+1:j+cols)
310P.rightCols<cols>() // P(:, end-cols+1:end)
311P.rightCols(cols) // P(:, end-cols+1:end)
312P.topRows<rows>() // P(1:rows, :)
313P.topRows(rows) // P(1:rows, :)
314P.middleRows<rows>(i) // P(i+1:i+rows, :)
315P.middleRows(i, rows) // P(i+1:i+rows, :)
316P.bottomRows<rows>() // P(end-rows+1:end, :)
317P.bottomRows(rows) // P(end-rows+1:end, :)
318P.topLeftCorner(rows, cols) // P(1:rows, 1:cols)
319P.topRightCorner(rows, cols) // P(1:rows, end-cols+1:end)
320P.bottomLeftCorner(rows, cols) // P(end-rows+1:end, 1:cols)
321P.bottomRightCorner(rows, cols) // P(end-rows+1:end, end-cols+1:end)
322P.topLeftCorner<rows,cols>() // P(1:rows, 1:cols)
323P.topRightCorner<rows,cols>() // P(1:rows, end-cols+1:end)
324P.bottomLeftCorner<rows,cols>() // P(end-rows+1:end, 1:cols)
325P.bottomRightCorner<rows,cols>() // P(end-rows+1:end, end-cols+1:end)
326
327
328
329
330// Of particular note is Eigen's swap function which is highly optimized.
331// Eigen // Matlab
332R.row(i) = P.col(j); // R(i, :) = P(:, i)
333R.col(j1).swap(mat1.col(j2)); // R(:, [j1 j2]) = R(:, [j2, j1])
334
335// Views, transpose, etc; all read-write except for .adjoint().
336// Eigen // Matlab
337R.adjoint() // R'
338R.transpose() // R.' or conj(R')
339R.diagonal() // diag(R)
340x.asDiagonal() // diag(x)
341R.transpose().colwise().reverse(); // rot90(R)
342R.conjugate() // conj(R)
343
344// All the same as Matlab, but matlab doesn't have *= style operators.
345// Matrix-vector. Matrix-matrix. Matrix-scalar.
346y = M*x; R = P*Q; R = P*s;
347a = b*M; R = P - Q; R = s*P;
348a *= M; R = P + Q; R = P/s;
349 R *= Q; R = s*P;
350 R += Q; R *= s;
351 R -= Q; R /= s;
352
353
354// Vectorized operations on each element independently
355// Eigen // Matlab
356R = P.cwiseProduct(Q); // R = P .* Q
357R = P.array() * s.array();// R = P .* s
358R = P.cwiseQuotient(Q); // R = P ./ Q
359R = P.array() / Q.array();// R = P ./ Q
360R = P.array() + s.array();// R = P + s
361R = P.array() - s.array();// R = P - s
362R.array() += s; // R = R + s
363R.array() -= s; // R = R - s
364R.array() < Q.array(); // R < Q
365R.array() <= Q.array(); // R <= Q
366R.cwiseInverse(); // 1 ./ P
367R.array().inverse(); // 1 ./ P
368R.array().sin() // sin(P)
369R.array().cos() // cos(P)
370R.array().pow(s) // P .^ s
371R.array().square() // P .^ 2
372R.array().cube() // P .^ 3
373R.cwiseSqrt() // sqrt(P)
374R.array().sqrt() // sqrt(P)
375R.array().exp() // exp(P)
376R.array().log() // log(P)
377R.cwiseMax(P) // max(R, P)
378R.array().max(P.array()) // max(R, P)
379R.cwiseMin(P) // min(R, P)
380R.array().min(P.array()) // min(R, P)
381R.cwiseAbs() //std::abs(P)
382R.array().abs() //std::abs(P)
383R.cwiseAbs2() //std::abs(P.^2)
384R.array().abs2() //std::abs(P.^2)
385(R.array() < s).select(P,Q); // (R < s ? P : Q)
386
387// Reductions.
388int r, c;
389// Eigen // Matlab
390R.minCoeff() // min(R(:))
391R.maxCoeff() // max(R(:))
392s = R.minCoeff(&r, &c) // [s, i] = min(R(:)); [r, c] = ind2sub(size(R), i);
393s = R.maxCoeff(&r, &c) // [s, i] = max(R(:)); [r, c] = ind2sub(size(R), i);
394R.sum() // sum(R(:))
395R.colwise().sum() // sum(R)
396R.rowwise().sum() // sum(R, 2) or sum(R')'
397R.prod() // prod(R(:))
398R.colwise().prod() // prod(R)
399R.rowwise().prod() // prod(R, 2) or prod(R')'
400R.trace() // trace(R)
401R.all() // all(R(:))
402R.colwise().all() // all(R)
403R.rowwise().all() // all(R, 2)
404R.any() // any(R(:))
405R.colwise().any() // any(R)
406R.rowwise().any() // any(R, 2)
407
408
409
410
411// Dot products, norms, etc.
412// Eigen // Matlab
413x.norm() // norm(x). Note that norm(R) doesn't work in Eigen.
414x.squaredNorm() // dot(x, x) Note the equivalence is not true for complex
415x.dot(y) // dot(x, y)
416x.cross(y) // cross(x, y) Requires #include <Eigen/Geometry>
417
419// Eigen // Matlab
420A.cast<kids_real>(); // kids_real(A)
421A.cast<float>(); // single(A)
422A.cast<int>(); // int32(A)
423A.real(); // real(A)
424A.imag(); // imag(A)
425// if the original type equals destination type, no work is done
426
427
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()
439
440// Eigenvalue problems
441// Eigen // Matlab
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<>
447
448
449// Use Map of Eigen, associate Map class with C++ pointer
450
451MapMd W(pW,1,N/F);
452std::cout << std::setiosflags(std::ios::scientific)
453 << std::setprecision(8) << W << std::endl;
454
455*****************************************/
456
457
458#endif // KIDS_LINALG_TPL_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
EigMX< kids_real > EigMXr
Definition linalg_tpl.h:55
EigMX< kids_complex > EigAXc
Definition linalg_tpl.h:58
Eigen::Map< EigMXr > MapMXr
Definition linalg_tpl.h:61
Eigen::Map< EigMXc > MapMXc
Definition linalg_tpl.h:62
Eigen::Map< EigAXr > MapAXr
Definition linalg_tpl.h:63
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, EigMajor > EigMX
Definition linalg_tpl.h:48
Eigen::Map< EigVXr > MapVXr
Definition linalg_tpl.h:59
EigMX< kids_real > EigAXr
Definition linalg_tpl.h:57
Eigen::Map< EigAXc > MapAXc
Definition linalg_tpl.h:64
EigVX< kids_complex > EigVXc
Definition linalg_tpl.h:54
Eigen::Matrix< T, Eigen::Dynamic, 1, EigMajor > EigVX
Definition linalg_tpl.h:45
Eigen::Array< T, Eigen::Dynamic, Eigen::Dynamic, EigMajor > EigAX
Definition linalg_tpl.h:51
EigVX< kids_real > EigVXr
Definition linalg_tpl.h:53
EigMX< kids_complex > EigMXc
Definition linalg_tpl.h:56
Eigen::Map< EigVXc > MapVXc
Definition linalg_tpl.h:60
< 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_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_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 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