MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
densemat.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_DENSEMAT
13#define MFEM_DENSEMAT
14
15#include "../config/config.hpp"
17#include "matrix.hpp"
18
19namespace mfem
20{
21
22/// Data type dense matrix using column-major storage
23class DenseMatrix : public Matrix
24{
25 friend class DenseTensor;
26 friend class DenseMatrixInverse;
27
28private:
29 Memory<real_t> data;
30
31 void Eigensystem(Vector &ev, DenseMatrix *evect = NULL);
32
33 void Eigensystem(DenseMatrix &b, Vector &ev, DenseMatrix *evect = NULL);
34
35 // Auxiliary method used in FNorm2() and FNorm()
36 void FNorm(real_t &scale_factor, real_t &scaled_fnorm2) const;
37
38public:
39 /** Default constructor for DenseMatrix.
40 Sets data = NULL and height = width = 0. */
42
43 /// Copy constructor
44 DenseMatrix(const DenseMatrix &);
45
46 /// Creates square matrix of size s.
47 explicit DenseMatrix(int s);
48
49 /// Creates rectangular matrix of size m x n.
50 DenseMatrix(int m, int n);
51
52 /// Creates rectangular matrix equal to the transpose of mat.
53 DenseMatrix(const DenseMatrix &mat, char ch);
54
55 /// Construct a DenseMatrix using an existing data array.
56 /** The DenseMatrix does not assume ownership of the data array, i.e. it will
57 not delete the array. */
58 DenseMatrix(real_t *d, int h, int w)
59 : Matrix(h, w) { UseExternalData(d, h, w); }
60
61 /// Create a dense matrix using a braced initializer list
62 /// The inner lists correspond to rows of the matrix
63 template <int M, int N, typename T = real_t>
64 explicit DenseMatrix(const T (&values)[M][N]) : DenseMatrix(M, N)
65 {
66 // DenseMatrix is column-major so copies have to be element-wise
67 for (int i = 0; i < M; i++)
68 {
69 for (int j = 0; j < N; j++)
70 {
71 (*this)(i,j) = values[i][j];
72 }
73 }
74 }
75
76 /// Change the data array and the size of the DenseMatrix.
77 /** The DenseMatrix does not assume ownership of the data array, i.e. it will
78 not delete the data array @a d. This method should not be used with
79 DenseMatrix that owns its current data array. */
80 void UseExternalData(real_t *d, int h, int w)
81 {
82 data.Wrap(d, h*w, false);
83 height = h; width = w;
84 }
85
86 /// Change the data array and the size of the DenseMatrix.
87 /** The DenseMatrix does not assume ownership of the data array, i.e. it will
88 not delete the new array @a d. This method will delete the current data
89 array, if owned. */
90 void Reset(real_t *d, int h, int w)
91 { if (OwnsData()) { data.Delete(); } UseExternalData(d, h, w); }
92
93 /** Clear the data array and the dimensions of the DenseMatrix. This method
94 should not be used with DenseMatrix that owns its current data array. */
95 void ClearExternalData() { data.Reset(); height = width = 0; }
96
97 /// Delete the matrix data array (if owned) and reset the matrix state.
98 void Clear()
99 { if (OwnsData()) { data.Delete(); } ClearExternalData(); }
100
101 /// For backward compatibility define Size to be synonym of Width()
102 int Size() const { return Width(); }
103
104 /// Change the size of the DenseMatrix to s x s.
105 void SetSize(int s) { SetSize(s, s); }
106
107 /// Change the size of the DenseMatrix to h x w.
108 void SetSize(int h, int w);
109
110 /// Returns the matrix data array.
111 inline real_t *Data() const
112 { return const_cast<real_t*>((const real_t*)data);}
113
114 /// Returns the matrix data array.
115 inline real_t *GetData() const { return Data(); }
116
117 Memory<real_t> &GetMemory() { return data; }
118 const Memory<real_t> &GetMemory() const { return data; }
119
120 /// Return the DenseMatrix data (host pointer) ownership flag.
121 inline bool OwnsData() const { return data.OwnsHostPtr(); }
122
123 /// Returns reference to a_{ij}.
124 inline real_t &operator()(int i, int j);
125
126 /// Returns constant reference to a_{ij}.
127 inline const real_t &operator()(int i, int j) const;
128
129 /// Matrix inner product: tr(A^t B)
130 real_t operator*(const DenseMatrix &m) const;
131
132 /// Trace of a square matrix
133 real_t Trace() const;
134
135 /// Returns reference to a_{ij}.
136 virtual real_t &Elem(int i, int j);
137
138 /// Returns constant reference to a_{ij}.
139 virtual const real_t &Elem(int i, int j) const;
140
141 /// Matrix vector multiplication.
142 void Mult(const real_t *x, real_t *y) const;
143
144 /// Matrix vector multiplication.
145 void Mult(const real_t *x, Vector &y) const;
146
147 /// Matrix vector multiplication.
148 void Mult(const Vector &x, real_t *y) const;
149
150 /// Matrix vector multiplication.
151 virtual void Mult(const Vector &x, Vector &y) const;
152
153 /// Multiply a vector with the transpose matrix.
154 void MultTranspose(const real_t *x, real_t *y) const;
155
156 /// Multiply a vector with the transpose matrix.
157 void MultTranspose(const real_t *x, Vector &y) const;
158
159 /// Multiply a vector with the transpose matrix.
160 void MultTranspose(const Vector &x, real_t *y) const;
161
162 /// Multiply a vector with the transpose matrix.
163 virtual void MultTranspose(const Vector &x, Vector &y) const;
164
165 using Operator::Mult;
167
168 /// y += a * A.x
169 virtual void AddMult(const Vector &x, Vector &y, const real_t a = 1.0) const;
170
171 /// y += a * A^t x
172 virtual void AddMultTranspose(const Vector &x, Vector &y,
173 const real_t a = 1.0) const;
174
175 /// y += a * A.x
176 void AddMult_a(real_t a, const Vector &x, Vector &y) const;
177
178 /// y += a * A^t x
179 void AddMultTranspose_a(real_t a, const Vector &x, Vector &y) const;
180
181 /// Compute y^t A x
182 real_t InnerProduct(const real_t *x, const real_t *y) const;
183
184 /// LeftScaling this = diag(s) * this
185 void LeftScaling(const Vector & s);
186 /// InvLeftScaling this = diag(1./s) * this
187 void InvLeftScaling(const Vector & s);
188 /// RightScaling: this = this * diag(s);
189 void RightScaling(const Vector & s);
190 /// InvRightScaling: this = this * diag(1./s);
191 void InvRightScaling(const Vector & s);
192 /// SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
193 void SymmetricScaling(const Vector & s);
194 /// InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
195 void InvSymmetricScaling(const Vector & s);
196
197 /// Compute y^t A x
198 real_t InnerProduct(const Vector &x, const Vector &y) const
199 { return InnerProduct(x.GetData(), y.GetData()); }
200
201 /// Returns a pointer to the inverse matrix.
202 virtual MatrixInverse *Inverse() const;
203
204 /// Replaces the current matrix with its inverse
205 void Invert();
206
207 /// Replaces the current matrix with its square root inverse
208 void SquareRootInverse();
209
210 /// Calculates the determinant of the matrix
211 /// (optimized for 2x2, 3x3, and 4x4 matrices)
212 real_t Det() const;
213
214 real_t Weight() const;
215
216 /** @brief Set the matrix to alpha * A, assuming that A has the same
217 dimensions as the matrix and uses column-major layout. */
218 void Set(real_t alpha, const real_t *A);
219 /// Set the matrix to alpha * A.
220 void Set(real_t alpha, const DenseMatrix &A)
221 {
222 SetSize(A.Height(), A.Width());
223 Set(alpha, A.GetData());
224 }
225
226 /// Adds the matrix A multiplied by the number c to the matrix.
227 void Add(const real_t c, const DenseMatrix &A);
228
229 /// Adds the matrix A multiplied by the number c to the matrix,
230 /// assuming A has the same dimensions and uses column-major layout.
231 void Add(const real_t c, const real_t *A);
232
233 /// Sets the matrix elements equal to constant c
235
236 /// Copy the matrix entries from the given array
237 DenseMatrix &operator=(const real_t *d);
238
239 /// Sets the matrix size and elements equal to those of m
241
242 DenseMatrix &operator+=(const real_t *m);
244
246
248
249 /// (*this) = -(*this)
250 void Neg();
251
252 /// Take the 2-norm of the columns of A and store in v
253 void Norm2(real_t *v) const;
254
255 /// Take the 2-norm of the columns of A and store in v
256 void Norm2(Vector &v) const
257 {
258 MFEM_ASSERT(v.Size() == Width(), "incompatible Vector size!");
259 Norm2(v.GetData());
260 }
261
262 /// Compute the norm ||A|| = max_{ij} |A_{ij}|
263 real_t MaxMaxNorm() const;
264
265 /// Compute the Frobenius norm of the matrix
266 real_t FNorm() const { real_t s, n2; FNorm(s, n2); return s*sqrt(n2); }
267
268 /// Compute the square of the Frobenius norm of the matrix
269 real_t FNorm2() const { real_t s, n2; FNorm(s, n2); return s*s*n2; }
270
271 /// Compute eigenvalues of A x = ev x where A = *this
273 { Eigensystem(ev); }
274
275 /// Compute eigenvalues and eigenvectors of A x = ev x where A = *this
277 { Eigensystem(ev, &evect); }
278
279 /// Compute eigenvalues and eigenvectors of A x = ev x where A = *this
281 { Eigensystem(ev, &evect); }
282
283 /** Compute generalized eigenvalues and eigenvectors of A x = ev B x,
284 where A = *this */
286 { Eigensystem(b, ev); }
287
288 /// Compute generalized eigenvalues of A x = ev B x, where A = *this
290 { Eigensystem(b, ev, &evect); }
291
292 /** Compute generalized eigenvalues and eigenvectors of A x = ev B x,
293 where A = *this */
295 { Eigensystem(b, ev, &evect); }
296
297 void SingularValues(Vector &sv) const;
298 int Rank(real_t tol) const;
299
300 /// Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
301 real_t CalcSingularvalue(const int i) const;
302
303 /** Return the eigenvalues (in increasing order) and eigenvectors of a
304 2x2 or 3x3 symmetric matrix. */
305 void CalcEigenvalues(real_t *lambda, real_t *vec) const;
306
307 void GetRow(int r, Vector &row) const;
308 void GetColumn(int c, Vector &col) const;
309 real_t *GetColumn(int col) { return data + col*height; }
310 const real_t *GetColumn(int col) const { return data + col*height; }
311
312 void GetColumnReference(int c, Vector &col)
313 { col.SetDataAndSize(data + c * height, height); }
314
315 void SetRow(int r, const real_t* row);
316 void SetRow(int r, const Vector &row);
317
318 void SetCol(int c, const real_t* col);
319 void SetCol(int c, const Vector &col);
320
321
322 /// Set all entries of a row to the specified value.
323 void SetRow(int row, real_t value);
324 /// Set all entries of a column to the specified value.
325 void SetCol(int col, real_t value);
326
327 /// Returns the diagonal of the matrix
328 void GetDiag(Vector &d) const;
329 /// Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|
330 void Getl1Diag(Vector &l) const;
331 /// Compute the row sums of the DenseMatrix
332 void GetRowSums(Vector &l) const;
333
334 /// Creates n x n diagonal matrix with diagonal elements c
335 void Diag(real_t c, int n);
336 /// Creates n x n diagonal matrix with diagonal given by diag
337 void Diag(real_t *diag, int n);
338
339 /// (*this) = (*this)^t
340 void Transpose();
341 /// (*this) = A^t
342 void Transpose(const DenseMatrix &A);
343 /// (*this) = 1/2 ((*this) + (*this)^t)
344 void Symmetrize();
345
346 void Lump();
347
348 /** Given a DShape matrix (from a scalar FE), stored in *this, returns the
349 CurlShape matrix. If *this is a N by D matrix, then curl is a D*N by
350 D*(D-1)/2 matrix. The size of curl must be set outside. The dimension D
351 can be either 2 or 3. In 2D this computes the scalar-valued curl of a
352 2D vector field */
353 void GradToCurl(DenseMatrix &curl);
354 /** Given a DShape matrix (from a scalar FE), stored in *this, returns the
355 CurlShape matrix. This computes the vector-valued curl of a scalar field.
356 *this is N by 2 matrix and curl is N by 2 matrix as well. */
358 /** Given a DShape matrix (from a scalar FE), stored in *this,
359 returns the DivShape vector. If *this is a N by dim matrix,
360 then div is a dim*N vector. The size of div must be set
361 outside. */
362 void GradToDiv(Vector &div);
363
364 /// Copy rows row1 through row2 from A to *this
365 void CopyRows(const DenseMatrix &A, int row1, int row2);
366 /// Copy columns col1 through col2 from A to *this
367 void CopyCols(const DenseMatrix &A, int col1, int col2);
368 /// Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this
369 void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco);
370 /// Copy matrix A to the location in *this at row_offset, col_offset
371 void CopyMN(const DenseMatrix &A, int row_offset, int col_offset);
372 /// Copy matrix A^t to the location in *this at row_offset, col_offset
373 void CopyMNt(const DenseMatrix &A, int row_offset, int col_offset);
374 /** Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this at
375 row_offset, col_offset */
376 void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco,
377 int row_offset, int col_offset);
378 /// Copy c on the diagonal of size n to *this at row_offset, col_offset
379 void CopyMNDiag(real_t c, int n, int row_offset, int col_offset);
380 /// Copy diag on the diagonal of size n to *this at row_offset, col_offset
381 void CopyMNDiag(real_t *diag, int n, int row_offset, int col_offset);
382 /// Copy All rows and columns except m and n from A
383 void CopyExceptMN(const DenseMatrix &A, int m, int n);
384
385 /// Perform (ro+i,co+j)+=A(i,j) for 0<=i<A.Height, 0<=j<A.Width
386 void AddMatrix(DenseMatrix &A, int ro, int co);
387 /// Perform (ro+i,co+j)+=a*A(i,j) for 0<=i<A.Height, 0<=j<A.Width
388 void AddMatrix(real_t a, const DenseMatrix &A, int ro, int co);
389
390 /** Get the square submatrix which corresponds to the given indices @a idx.
391 Note: the @a A matrix will be resized to accommodate the data */
392 void GetSubMatrix(const Array<int> & idx, DenseMatrix & A) const;
393
394 /** Get the rectangular submatrix which corresponds to the given indices
395 @a idx_i and @a idx_j. Note: the @a A matrix will be resized to
396 accommodate the data */
397 void GetSubMatrix(const Array<int> & idx_i, const Array<int> & idx_j,
398 DenseMatrix & A) const;
399
400 /** Get the square submatrix which corresponds to the range
401 [ @a ibeg, @a iend ). Note: the @a A matrix will be resized
402 to accommodate the data */
403 void GetSubMatrix(int ibeg, int iend, DenseMatrix & A);
404
405 /** Get the square submatrix which corresponds to the range
406 i ∈ [ @a ibeg, @a iend ) and j ∈ [ @a jbeg, @a jend )
407 Note: the @a A matrix will be resized to accommodate the data */
408 void GetSubMatrix(int ibeg, int iend, int jbeg, int jend, DenseMatrix & A);
409
410 /// Set (*this)(idx[i],idx[j]) = A(i,j)
411 void SetSubMatrix(const Array<int> & idx, const DenseMatrix & A);
412
413 /// Set (*this)(idx_i[i],idx_j[j]) = A(i,j)
414 void SetSubMatrix(const Array<int> & idx_i, const Array<int> & idx_j,
415 const DenseMatrix & A);
416
417 /** Set a submatrix of (this) to the given matrix @a A
418 with row and column offset @a ibeg */
419 void SetSubMatrix(int ibeg, const DenseMatrix & A);
420
421 /** Set a submatrix of (this) to the given matrix @a A
422 with row and column offset @a ibeg and @a jbeg respectively */
423 void SetSubMatrix(int ibeg, int jbeg, const DenseMatrix & A);
424
425 /// (*this)(idx[i],idx[j]) += A(i,j)
426 void AddSubMatrix(const Array<int> & idx, const DenseMatrix & A);
427
428 /// (*this)(idx_i[i],idx_j[j]) += A(i,j)
429 void AddSubMatrix(const Array<int> & idx_i, const Array<int> & idx_j,
430 const DenseMatrix & A);
431
432 /** Add the submatrix @a A to this with row and column offset @a ibeg */
433 void AddSubMatrix(int ibeg, const DenseMatrix & A);
434
435 /** Add the submatrix @a A to this with row and column offsets
436 @a ibeg and @a jbeg respectively */
437 void AddSubMatrix(int ibeg, int jbeg, const DenseMatrix & A);
438
439 /// Add the matrix 'data' to the Vector 'v' at the given 'offset'
440 void AddToVector(int offset, Vector &v) const;
441 /// Get the matrix 'data' from the Vector 'v' at the given 'offset'
442 void GetFromVector(int offset, const Vector &v);
443 /** If (dofs[i] < 0 and dofs[j] >= 0) or (dofs[i] >= 0 and dofs[j] < 0)
444 then (*this)(i,j) = -(*this)(i,j). */
445 void AdjustDofDirection(Array<int> &dofs);
446
447 /// Replace small entries, abs(a_ij) <= eps, with zero.
448 void Threshold(real_t eps);
449
450 /** Count the number of entries in the matrix for which isfinite
451 is false, i.e. the entry is a NaN or +/-Inf. */
453
454 /// Prints matrix to stream out.
455 virtual void Print(std::ostream &out = mfem::out, int width_ = 4) const;
456 virtual void PrintMatlab(std::ostream &out = mfem::out) const;
457 /// Prints the transpose matrix to stream out.
458 virtual void PrintT(std::ostream &out = mfem::out, int width_ = 4) const;
459
460 /// Invert and print the numerical conditioning of the inversion.
461 void TestInversion();
462
463 std::size_t MemoryUsage() const { return data.Capacity() * sizeof(real_t); }
464
465 /// Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
466 const real_t *Read(bool on_dev = true) const
467 { return mfem::Read(data, Height()*Width(), on_dev); }
468
469 /// Shortcut for mfem::Read(GetMemory(), TotalSize(), false).
470 const real_t *HostRead() const
471 { return mfem::Read(data, Height()*Width(), false); }
472
473 /// Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
474 real_t *Write(bool on_dev = true)
475 { return mfem::Write(data, Height()*Width(), on_dev); }
476
477 /// Shortcut for mfem::Write(GetMemory(), TotalSize(), false).
479 { return mfem::Write(data, Height()*Width(), false); }
480
481 /// Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev).
482 real_t *ReadWrite(bool on_dev = true)
483 { return mfem::ReadWrite(data, Height()*Width(), on_dev); }
484
485 /// Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false).
487 { return mfem::ReadWrite(data, Height()*Width(), false); }
488
489 void Swap(DenseMatrix &other);
490
491 /// Destroys dense matrix.
492 virtual ~DenseMatrix();
493};
494
495/// C = A + alpha*B
496void Add(const DenseMatrix &A, const DenseMatrix &B,
497 real_t alpha, DenseMatrix &C);
498
499/// C = alpha*A + beta*B
500void Add(real_t alpha, const real_t *A,
501 real_t beta, const real_t *B, DenseMatrix &C);
502
503/// C = alpha*A + beta*B
504void Add(real_t alpha, const DenseMatrix &A,
505 real_t beta, const DenseMatrix &B, DenseMatrix &C);
506
507/// @brief Solves the dense linear system, `A * X = B` for `X`
508///
509/// @param [in,out] A the square matrix for the linear system
510/// @param [in,out] X the rhs vector, B, on input, the solution, X, on output.
511/// @param [in] TOL optional fuzzy comparison tolerance. Defaults to 1e-9.
512///
513/// @return status set to true if successful, otherwise, false.
514///
515/// @note This routine may replace the contents of the input Matrix, A, with the
516/// corresponding LU factorization of the matrix. Matrices of size 1x1 and
517/// 2x2 are handled explicitly.
518///
519/// @pre A.IsSquare() == true
520/// @pre X != nullptr
521bool LinearSolve(DenseMatrix& A, real_t* X, real_t TOL = 1.e-9);
522
523/// Matrix matrix multiplication. A = B * C.
524void Mult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a);
525
526/// Matrix matrix multiplication. A += B * C.
527void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a);
528
529/// Matrix matrix multiplication. A += alpha * B * C.
530void AddMult_a(real_t alpha, const DenseMatrix &b, const DenseMatrix &c,
531 DenseMatrix &a);
532
533/** Calculate the adjugate of a matrix (for NxN matrices, N=1,2,3) or the matrix
534 adj(A^t.A).A^t for rectangular matrices (2x1, 3x1, or 3x2). This operation
535 is well defined even when the matrix is not full rank. */
536void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja);
537
538/// Calculate the transposed adjugate of a matrix (for NxN matrices, N=1,2,3)
539void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat);
540
541/** Calculate the inverse of a matrix (for NxN matrices, N=1,2,3) or the
542 left inverse (A^t.A)^{-1}.A^t (for 2x1, 3x1, or 3x2 matrices) */
543void CalcInverse(const DenseMatrix &a, DenseMatrix &inva);
544
545/// Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
546void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva);
547
548/** For a given Nx(N-1) (N=2,3) matrix J, compute a vector n such that
549 n_k = (-1)^{k+1} det(J_k), k=1,..,N, where J_k is the matrix J with the
550 k-th row removed. Note: J^t.n = 0, det([n|J])=|n|^2=det(J^t.J). */
551void CalcOrtho(const DenseMatrix &J, Vector &n);
552
553/// Calculate the matrix A.At
554void MultAAt(const DenseMatrix &a, DenseMatrix &aat);
555
556/// ADAt = A D A^t, where D is diagonal
557void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt);
558
559/// ADAt += A D A^t, where D is diagonal
560void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt);
561
562/// Multiply a matrix A with the transpose of a matrix B: A*Bt
563void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt);
564
565/// ADBt = A D B^t, where D is diagonal
566void MultADBt(const DenseMatrix &A, const Vector &D,
567 const DenseMatrix &B, DenseMatrix &ADBt);
568
569/// ABt += A * B^t
570void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt);
571
572/// ADBt = A D B^t, where D is diagonal
573void AddMultADBt(const DenseMatrix &A, const Vector &D,
574 const DenseMatrix &B, DenseMatrix &ADBt);
575
576/// ABt += a * A * B^t
577void AddMult_a_ABt(real_t a, const DenseMatrix &A, const DenseMatrix &B,
578 DenseMatrix &ABt);
579
580/// Multiply the transpose of a matrix A with a matrix B: At*B
581void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB);
582
583/// AtB += A^t * B
584void AddMultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB);
585
586/// AtB += a * A^t * B
587void AddMult_a_AtB(real_t a, const DenseMatrix &A, const DenseMatrix &B,
588 DenseMatrix &AtB);
589
590/// AAt += a * A * A^t
591void AddMult_a_AAt(real_t a, const DenseMatrix &A, DenseMatrix &AAt);
592
593/// AAt = a * A * A^t
594void Mult_a_AAt(real_t a, const DenseMatrix &A, DenseMatrix &AAt);
595
596/// Make a matrix from a vector V.Vt
597void MultVVt(const Vector &v, DenseMatrix &vvt);
598
599void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt);
600
601/// VWt += v w^t
602void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt);
603
604/// VVt += v v^t
605void AddMultVVt(const Vector &v, DenseMatrix &VWt);
606
607/// VWt += a * v w^t
608void AddMult_a_VWt(const real_t a, const Vector &v, const Vector &w,
609 DenseMatrix &VWt);
610
611/// VVt += a * v v^t
612void AddMult_a_VVt(const real_t a, const Vector &v, DenseMatrix &VVt);
613
614/** Computes matrix P^t * A * P. Note: The @a RAP matrix will be resized
615 to accommodate the data */
616void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix & RAP);
617
618/** Computes the matrix Rt^t * A * P. Note: The @a RAP matrix will be resized
619 to accommodate the data */
620void RAP(const DenseMatrix &Rt, const DenseMatrix &A,
621 const DenseMatrix &P, DenseMatrix & RAP);
622
623/** Abstract class that can compute factorization of external data and perform various
624 operations with the factored data. */
626{
627public:
628
630
632
633 Factors(real_t *data_) : data(data_) { }
634
635 virtual bool Factor(int m, real_t TOL = 0.0)
636 {
637 mfem_error("Factors::Factors(...)");
638 return false;
639 }
640
641 virtual real_t Det(int m) const
642 {
643 mfem_error("Factors::Det(...)");
644 return 0.;
645 }
646
647 virtual void Solve(int m, int n, real_t *X) const
648 {
649 mfem_error("Factors::Solve(...)");
650 }
651
652 virtual void GetInverseMatrix(int m, real_t *X) const
653 {
654 mfem_error("Factors::GetInverseMatrix(...)");
655 }
656
657 virtual ~Factors() {}
658};
659
660
661/** Class that can compute LU factorization of external data and perform various
662 operations with the factored data. */
663class LUFactors : public Factors
664{
665public:
666 int *ipiv;
667#ifdef MFEM_USE_LAPACK
668 static const int ipiv_base = 1;
669#else
670 static const int ipiv_base = 0;
671#endif
672
673 /** With this constructor, the (public) data and ipiv members should be set
674 explicitly before calling class methods. */
676
677 LUFactors(real_t *data_, int *ipiv_) : Factors(data_), ipiv(ipiv_) { }
678
679 /**
680 * @brief Compute the LU factorization of the current matrix
681 *
682 * Factorize the current matrix of size (m x m) overwriting it with the
683 * LU factors. The factorization is such that L.U = P.A, where A is the
684 * original matrix and P is a permutation matrix represented by ipiv.
685 *
686 * @param [in] m size of the square matrix
687 * @param [in] TOL optional fuzzy comparison tolerance. Defaults to 0.0.
688 *
689 * @return status set to true if successful, otherwise, false.
690 */
691 virtual bool Factor(int m, real_t TOL = 0.0);
692
693 /** Assuming L.U = P.A factored data of size (m x m), compute |A|
694 from the diagonal values of U and the permutation information. */
695 virtual real_t Det(int m) const;
696
697 /** Assuming L.U = P.A factored data of size (m x m), compute X <- A X,
698 for a matrix X of size (m x n). */
699 void Mult(int m, int n, real_t *X) const;
700
701 /** Assuming L.U = P.A factored data of size (m x m), compute
702 X <- L^{-1} P X, for a matrix X of size (m x n). */
703 void LSolve(int m, int n, real_t *X) const;
704
705 /** Assuming L.U = P.A factored data of size (m x m), compute
706 X <- U^{-1} X, for a matrix X of size (m x n). */
707 void USolve(int m, int n, real_t *X) const;
708
709 /** Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1} X,
710 for a matrix X of size (m x n). */
711 virtual void Solve(int m, int n, real_t *X) const;
712
713 /** Assuming L.U = P.A factored data of size (m x m), compute X <- X A^{-1},
714 for a matrix X of size (n x m). */
715 void RightSolve(int m, int n, real_t *X) const;
716
717 /// Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1}.
718 virtual void GetInverseMatrix(int m, real_t *X) const;
719
720 /** Given an (n x m) matrix A21, compute X2 <- X2 - A21 X1, for matrices X1,
721 and X2 of size (m x r) and (n x r), respectively. */
722 static void SubMult(int m, int n, int r, const real_t *A21,
723 const real_t *X1, real_t *X2);
724
725 /** Assuming P.A = L.U factored data of size (m x m), compute the 2x2 block
726 decomposition:
727 | P 0 | | A A12 | = | L 0 | | U U12 |
728 | 0 I | | A21 A22 | | L21 I | | 0 S22 |
729 where A12, A21, and A22 are matrices of size (m x n), (n x m), and
730 (n x n), respectively. The blocks are overwritten as follows:
731 A12 <- U12 = L^{-1} P A12
732 A21 <- L21 = A21 U^{-1}
733 A22 <- S22 = A22 - L21 U12.
734 The block S22 is the Schur complement. */
735 void BlockFactor(int m, int n, real_t *A12, real_t *A21, real_t *A22) const;
736
737 /** Given BlockFactor()'d data, perform the forward block solve for the
738 linear system:
739 | A A12 | | X1 | = | B1 |
740 | A21 A22 | | X2 | | B2 |
741 written in the factored form:
742 | L 0 | | U U12 | | X1 | = | P 0 | | B1 |
743 | L21 I | | 0 S22 | | X2 | | 0 I | | B2 |.
744 The resulting blocks Y1, Y2 solve the system:
745 | L 0 | | Y1 | = | P 0 | | B1 |
746 | L21 I | | Y2 | | 0 I | | B2 |
747 The blocks are overwritten as follows:
748 B1 <- Y1 = L^{-1} P B1
749 B2 <- Y2 = B2 - L21 Y1 = B2 - A21 A^{-1} B1
750 The blocks B1/Y1 and B2/Y2 are of size (m x r) and (n x r), respectively.
751 The Schur complement system is given by: S22 X2 = Y2. */
752 void BlockForwSolve(int m, int n, int r, const real_t *L21,
753 real_t *B1, real_t *B2) const;
754
755 /** Given BlockFactor()'d data, perform the backward block solve in
756 | U U12 | | X1 | = | Y1 |
757 | 0 S22 | | X2 | | Y2 |.
758 The input is the solution block X2 and the block Y1 resulting from
759 BlockForwSolve(). The result block X1 overwrites input block Y1:
760 Y1 <- X1 = U^{-1} (Y1 - U12 X2). */
761 void BlockBackSolve(int m, int n, int r, const real_t *U12,
762 const real_t *X2, real_t *Y1) const;
763};
764
765
766/** Class that can compute Cholesky factorizations of external data of an
767 SPD matrix and perform various operations with the factored data. */
769{
770public:
771
772 /** With this constructor, the (public) data should be set
773 explicitly before calling class methods. */
775
776 CholeskyFactors(real_t *data_) : Factors(data_) { }
777
778 /**
779 * @brief Compute the Cholesky factorization of the current matrix
780 *
781 * Factorize the current matrix of size (m x m) overwriting it with the
782 * Cholesky factors. The factorization is such that LL^t = A, where A is the
783 * original matrix
784 *
785 * @param [in] m size of the square matrix
786 * @param [in] TOL optional fuzzy comparison tolerance. Defaults to 0.0.
787 *
788 * @return status set to true if successful, otherwise, false.
789 */
790 virtual bool Factor(int m, real_t TOL = 0.0);
791
792 /** Assuming LL^t = A factored data of size (m x m), compute |A|
793 from the diagonal values of L */
794 virtual real_t Det(int m) const;
795
796 /** Assuming L.L^t = A factored data of size (m x m), compute X <- L X,
797 for a matrix X of size (m x n). */
798 void LMult(int m, int n, real_t *X) const;
799
800 /** Assuming L.L^t = A factored data of size (m x m), compute X <- L^t X,
801 for a matrix X of size (m x n). */
802 void UMult(int m, int n, real_t *X) const;
803
804 /** Assuming L L^t = A factored data of size (m x m), compute
805 X <- L^{-1} X, for a matrix X of size (m x n). */
806 void LSolve(int m, int n, real_t *X) const;
807
808 /** Assuming L L^t = A factored data of size (m x m), compute
809 X <- L^{-t} X, for a matrix X of size (m x n). */
810 void USolve(int m, int n, real_t *X) const;
811
812 /** Assuming L.L^t = A factored data of size (m x m), compute X <- A^{-1} X,
813 for a matrix X of size (m x n). */
814 virtual void Solve(int m, int n, real_t *X) const;
815
816 /** Assuming L.L^t = A factored data of size (m x m), compute X <- X A^{-1},
817 for a matrix X of size (n x m). */
818 void RightSolve(int m, int n, real_t *X) const;
819
820 /// Assuming L.L^t = A factored data of size (m x m), compute X <- A^{-1}.
821 virtual void GetInverseMatrix(int m, real_t *X) const;
822
823};
824
825
826/** Data type for inverse of square dense matrix.
827 Stores matrix factors, i.e., Cholesky factors if the matrix is SPD,
828 LU otherwise. */
830{
831private:
832 const DenseMatrix *a;
833 Factors * factors = nullptr;
834 bool spd = false;
835
836 void Init(int m);
837 bool own_data = false;
838public:
839 /// Default constructor.
840 DenseMatrixInverse(bool spd_=false) : a(NULL), spd(spd_) { Init(0); }
841
842 /** Creates square dense matrix. Computes factorization of mat
843 and stores its factors. */
844 DenseMatrixInverse(const DenseMatrix &mat, bool spd_ = false);
845
846 /// Same as above but does not factorize the matrix.
847 DenseMatrixInverse(const DenseMatrix *mat, bool spd_ = false);
848
849 /// Get the size of the inverse matrix
850 int Size() const { return Width(); }
851
852 /// Factor the current DenseMatrix, *a
853 void Factor();
854
855 /// Factor a new DenseMatrix of the same size
856 void Factor(const DenseMatrix &mat);
857
858 virtual void SetOperator(const Operator &op);
859
860 /// Matrix vector multiplication with the inverse of dense matrix.
861 void Mult(const real_t *x, real_t *y) const;
862
863 /// Matrix vector multiplication with the inverse of dense matrix.
864 virtual void Mult(const Vector &x, Vector &y) const;
865
866 /// Multiply the inverse matrix by another matrix: X = A^{-1} B.
867 void Mult(const DenseMatrix &B, DenseMatrix &X) const;
868
869 /// Multiply the inverse matrix by another matrix: X <- A^{-1} X.
870 void Mult(DenseMatrix &X) const {factors->Solve(width, X.Width(), X.Data());}
871
872 using Operator::Mult;
873
874 /// Compute and return the inverse matrix in Ainv.
875 void GetInverseMatrix(DenseMatrix &Ainv) const;
876
877 /// Compute the determinant of the original DenseMatrix using the LU factors.
878 real_t Det() const { return factors->Det(width); }
879
880 /// Print the numerical conditioning of the inversion: ||A^{-1} A - I||.
881 void TestInversion();
882
883 /// Destroys dense inverse matrix.
884 virtual ~DenseMatrixInverse();
885};
886
887#ifdef MFEM_USE_LAPACK
888
890{
891 DenseMatrix &mat;
892 Vector EVal;
893 DenseMatrix EVect;
894 Vector ev;
895 int n;
896 real_t *work;
897 char jobz, uplo;
898 int lwork, info;
899public:
900
903 void Eval();
904 Vector &Eigenvalues() { return EVal; }
905 DenseMatrix &Eigenvectors() { return EVect; }
906 real_t Eigenvalue(int i) { return EVal(i); }
907 const Vector &Eigenvector(int i)
908 {
909 ev.SetData(EVect.Data() + i * EVect.Height());
910 return ev;
911 }
913};
914
916{
917 DenseMatrix &A;
918 DenseMatrix &B;
919 DenseMatrix A_copy;
920 DenseMatrix B_copy;
921 Vector evalues_r;
922 Vector evalues_i;
923 DenseMatrix Vr;
924 DenseMatrix Vl;
925 int n;
926
927 real_t *alphar;
928 real_t *alphai;
929 real_t *beta;
930 real_t *work;
931 char jobvl, jobvr;
932 int lwork, info;
933
934public:
935
937 bool left_eigen_vectors = false,
938 bool right_eigen_vectors = false);
939 void Eval();
940 Vector &EigenvaluesRealPart() { return evalues_r; }
941 Vector &EigenvaluesImagPart() { return evalues_i; }
942 real_t EigenvalueRealPart(int i) { return evalues_r(i); }
943 real_t EigenvalueImagPart(int i) { return evalues_i(i); }
947};
948
949/**
950 @brief Class for Singular Value Decomposition of a DenseMatrix
951
952 Singular Value Decomposition (SVD) of a DenseMatrix with the use of the DGESVD
953 driver from LAPACK.
954 */
956{
957 DenseMatrix Mc;
958 Vector sv;
959 DenseMatrix U,Vt;
960 int m, n;
961
962#ifdef MFEM_USE_LAPACK
963 real_t *work;
964 char jobu, jobvt;
965 int lwork, info;
966#endif
967
968 void Init();
969public:
970
971 /**
972 @brief Constructor for the DenseMatrixSVD
973
974 Constructor for the DenseMatrixSVD with LAPACK. The parameters for the left
975 and right singular vectors can be choosen according to the parameters for
976 the LAPACK DGESVD.
977
978 @param [in] M matrix to set the size to n=M.Height(), m=M.Width()
979 @param [in] left_singular_vectors optional parameter to define if first
980 left singular vectors should be computed
981 @param [in] right_singular_vectors optional parameter to define if first
982 right singular vectors should be computed
983 */
984 MFEM_DEPRECATED DenseMatrixSVD(DenseMatrix &M,
985 bool left_singular_vectors=false,
986 bool right_singular_vectors=false);
987
988 /**
989 @brief Constructor for the DenseMatrixSVD
990
991 Constructor for the DenseMatrixSVD with LAPACK. The parameters for the left
992 and right singular
993 vectors can be choosen according to the parameters for the LAPACK DGESVD.
994
995 @param [in] h height of the matrix
996 @param [in] w width of the matrix
997 @param [in] left_singular_vectors optional parameter to define if first
998 left singular vectors should be computed
999 @param [in] right_singular_vectors optional parameter to define if first
1000 right singular vectors should be computed
1001 */
1002 MFEM_DEPRECATED DenseMatrixSVD(int h, int w,
1003 bool left_singular_vectors=false,
1004 bool right_singular_vectors=false);
1005
1006 /**
1007 @brief Constructor for the DenseMatrixSVD
1008
1009 Constructor for the DenseMatrixSVD with LAPACK. The parameters for the left
1010 and right singular vectors can be choosen according to the parameters for
1011 the LAPACK DGESVD.
1012
1013 @param [in] M matrix to set the size to n=M.Height(), m=M.Width()
1014 @param [in] left_singular_vectors optional parameter to define which left
1015 singular vectors should be computed
1016 @param [in] right_singular_vectors optional parameter to define which right
1017 singular vectors should be computed
1018
1019 Options for computation of singular vectors:
1020
1021 'A': All singular vectors are computed (default)
1022
1023 'S': The first min(n,m) singular vectors are computed
1024
1025 'N': No singular vectors are computed
1026 */
1028 char left_singular_vectors='A',
1029 char right_singular_vectors='A');
1030
1031 /**
1032 @brief Constructor for the DenseMatrixSVD
1033
1034 Constructor for the DenseMatrixSVD with LAPACK. The parameters for the left
1035 and right singular vectors can be choosen according to the
1036 parameters for the LAPACK DGESVD.
1037
1038 @param [in] h height of the matrix
1039 @param [in] w width of the matrix
1040 @param [in] left_singular_vectors optional parameter to define which left
1041 singular vectors should be computed
1042 @param [in] right_singular_vectors optional parameter to define which right
1043 singular vectors should be computed
1044
1045 Options for computation of singular vectors:
1046
1047 'A': All singular vectors are computed (default)
1048
1049 'S': The first min(n,m) singular vectors are computed
1050
1051 'N': No singular vectors are computed
1052 */
1053 DenseMatrixSVD(int h, int w,
1054 char left_singular_vectors='A',
1055 char right_singular_vectors='A');
1056
1057 /**
1058 @brief Evaluate the SVD
1059
1060 Call of the DGESVD driver from LAPACK for the DenseMatrix M. The singular
1061 vectors are computed according to the setup in the call of the constructor.
1062
1063 @param [in] M DenseMatrix the SVD should be evaluated for
1064 */
1065 void Eval(DenseMatrix &M);
1066
1067 /**
1068 @brief Return singular values
1069
1070 @return sv Vector containing all singular values
1071 */
1072 Vector &Singularvalues() { return sv; }
1073
1074 /**
1075 @brief Return specific singular value
1076
1077 @return sv(i) i-th singular value
1078 */
1079 real_t Singularvalue(int i) { return sv(i); }
1080
1081 /**
1082 @brief Return left singular vectors
1083
1084 @return U DenseMatrix containing left singular vectors
1085 */
1087
1088 /**
1089 @brief Return right singular vectors
1090
1091 @return Vt DenseMatrix containing right singular vectors
1092 */
1095};
1096
1097#endif // if MFEM_USE_LAPACK
1098
1099
1100class Table;
1101
1102/// Rank 3 tensor (array of matrices)
1104{
1105private:
1106 mutable DenseMatrix Mk;
1107 Memory<real_t> tdata;
1108 int nk;
1109
1110public:
1112 {
1113 nk = 0;
1114 }
1115
1116 DenseTensor(int i, int j, int k)
1117 : Mk(NULL, i, j)
1118 {
1119 nk = k;
1120 tdata.New(i*j*k);
1121 }
1122
1123 DenseTensor(real_t *d, int i, int j, int k)
1124 : Mk(NULL, i, j)
1125 {
1126 nk = k;
1127 tdata.Wrap(d, i*j*k, false);
1128 }
1129
1130 DenseTensor(int i, int j, int k, MemoryType mt)
1131 : Mk(NULL, i, j)
1132 {
1133 nk = k;
1134 tdata.New(i*j*k, mt);
1135 }
1136
1137 /// Copy constructor: deep copy
1139 : Mk(NULL, other.Mk.height, other.Mk.width), nk(other.nk)
1140 {
1141 const int size = Mk.Height()*Mk.Width()*nk;
1142 if (size > 0)
1143 {
1144 tdata.New(size, other.tdata.GetMemoryType());
1145 tdata.CopyFrom(other.tdata, size);
1146 }
1147 }
1148
1149 int SizeI() const { return Mk.Height(); }
1150 int SizeJ() const { return Mk.Width(); }
1151 int SizeK() const { return nk; }
1152
1153 int TotalSize() const { return SizeI()*SizeJ()*SizeK(); }
1154
1155 void SetSize(int i, int j, int k, MemoryType mt_ = MemoryType::PRESERVE)
1156 {
1157 const MemoryType mt = mt_ == MemoryType::PRESERVE ? tdata.GetMemoryType() : mt_;
1158 tdata.Delete();
1159 Mk.UseExternalData(NULL, i, j);
1160 nk = k;
1161 tdata.New(i*j*k, mt);
1162 }
1163
1164 void UseExternalData(real_t *ext_data, int i, int j, int k)
1165 {
1166 tdata.Delete();
1167 Mk.UseExternalData(NULL, i, j);
1168 nk = k;
1169 tdata.Wrap(ext_data, i*j*k, false);
1170 }
1171
1172 /// Sets the tensor elements equal to constant c
1174
1175 /// Copy assignment operator (performs a deep copy)
1176 DenseTensor &operator=(const DenseTensor &other);
1177
1179 {
1180 MFEM_ASSERT_INDEX_IN_RANGE(k, 0, SizeK());
1181 Mk.data = Memory<real_t>(GetData(k), SizeI()*SizeJ(), false);
1182 return Mk;
1183 }
1184 const DenseMatrix &operator()(int k) const
1185 {
1186 MFEM_ASSERT_INDEX_IN_RANGE(k, 0, SizeK());
1187 Mk.data = Memory<real_t>(const_cast<real_t*>(GetData(k)), SizeI()*SizeJ(),
1188 false);
1189 return Mk;
1190 }
1191
1192 real_t &operator()(int i, int j, int k)
1193 {
1194 MFEM_ASSERT_INDEX_IN_RANGE(i, 0, SizeI());
1195 MFEM_ASSERT_INDEX_IN_RANGE(j, 0, SizeJ());
1196 MFEM_ASSERT_INDEX_IN_RANGE(k, 0, SizeK());
1197 return tdata[i+SizeI()*(j+SizeJ()*k)];
1198 }
1199
1200 const real_t &operator()(int i, int j, int k) const
1201 {
1202 MFEM_ASSERT_INDEX_IN_RANGE(i, 0, SizeI());
1203 MFEM_ASSERT_INDEX_IN_RANGE(j, 0, SizeJ());
1204 MFEM_ASSERT_INDEX_IN_RANGE(k, 0, SizeK());
1205 return tdata[i+SizeI()*(j+SizeJ()*k)];
1206 }
1207
1209 {
1210 MFEM_ASSERT_INDEX_IN_RANGE(k, 0, SizeK());
1211 return tdata+k*Mk.Height()*Mk.Width();
1212 }
1213
1214 const real_t *GetData(int k) const
1215 {
1216 MFEM_ASSERT_INDEX_IN_RANGE(k, 0, SizeK());
1217 return tdata+k*Mk.Height()*Mk.Width();
1218 }
1219
1220 real_t *Data() { return tdata; }
1221
1222 const real_t *Data() const { return tdata; }
1223
1224 Memory<real_t> &GetMemory() { return tdata; }
1225 const Memory<real_t> &GetMemory() const { return tdata; }
1226
1227 /** Matrix-vector product from unassembled element matrices, assuming both
1228 'x' and 'y' use the same elem_dof table. */
1229 void AddMult(const Table &elem_dof, const Vector &x, Vector &y) const;
1230
1231 void Clear()
1232 { UseExternalData(NULL, 0, 0, 0); }
1233
1234 std::size_t MemoryUsage() const { return nk*Mk.MemoryUsage(); }
1235
1236 /// Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
1237 const real_t *Read(bool on_dev = true) const
1238 { return mfem::Read(tdata, Mk.Height()*Mk.Width()*nk, on_dev); }
1239
1240 /// Shortcut for mfem::Read(GetMemory(), TotalSize(), false).
1241 const real_t *HostRead() const
1242 { return mfem::Read(tdata, Mk.Height()*Mk.Width()*nk, false); }
1243
1244 /// Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
1245 real_t *Write(bool on_dev = true)
1246 { return mfem::Write(tdata, Mk.Height()*Mk.Width()*nk, on_dev); }
1247
1248 /// Shortcut for mfem::Write(GetMemory(), TotalSize(), false).
1250 { return mfem::Write(tdata, Mk.Height()*Mk.Width()*nk, false); }
1251
1252 /// Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev).
1253 real_t *ReadWrite(bool on_dev = true)
1254 { return mfem::ReadWrite(tdata, Mk.Height()*Mk.Width()*nk, on_dev); }
1255
1256 /// Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false).
1258 { return mfem::ReadWrite(tdata, Mk.Height()*Mk.Width()*nk, false); }
1259
1261 {
1262 mfem::Swap(tdata, t.tdata);
1263 mfem::Swap(nk, t.nk);
1264 Mk.Swap(t.Mk);
1265 }
1266
1267 ~DenseTensor() { tdata.Delete(); }
1268};
1269
1270/** @brief Compute the LU factorization of a batch of matrices
1271
1272 Factorize n matrices of size (m x m) stored in a dense tensor overwriting it
1273 with the LU factors. The factorization is such that L.U = Piv.A, where A is
1274 the original matrix and Piv is a permutation matrix represented by P.
1275
1276 @param [in, out] Mlu batch of square matrices - dimension m x m x n.
1277 @param [out] P array storing pivot information - dimension m x n.
1278 @param [in] TOL optional fuzzy comparison tolerance. Defaults to 0.0. */
1279void BatchLUFactor(DenseTensor &Mlu, Array<int> &P, const real_t TOL = 0.0);
1280
1281/** @brief Solve batch linear systems
1282
1283 Assuming L.U = P.A for n factored matrices (m x m), compute x <- A x, for n
1284 companion vectors.
1285
1286 @param [in] Mlu batch of LU factors for matrix M - dimension m x m x n.
1287 @param [in] P array storing pivot information - dimension m x n.
1288 @param [in, out] X vector storing right-hand side and then solution -
1289 dimension m x n. */
1290void BatchLUSolve(const DenseTensor &Mlu, const Array<int> &P, Vector &X);
1291
1292
1293// Inline methods
1294
1295inline real_t &DenseMatrix::operator()(int i, int j)
1296{
1297 MFEM_ASSERT(data && i >= 0 && i < height && j >= 0 && j < width, "");
1298 return data[i+j*height];
1299}
1300
1301inline const real_t &DenseMatrix::operator()(int i, int j) const
1302{
1303 MFEM_ASSERT(data && i >= 0 && i < height && j >= 0 && j < width, "");
1304 return data[i+j*height];
1305}
1306
1307} // namespace mfem
1308
1309#endif
void UMult(int m, int n, real_t *X) const
void RightSolve(int m, int n, real_t *X) const
virtual void GetInverseMatrix(int m, real_t *X) const
Assuming L.L^t = A factored data of size (m x m), compute X <- A^{-1}.
CholeskyFactors(real_t *data_)
Definition densemat.hpp:776
void USolve(int m, int n, real_t *X) const
virtual real_t Det(int m) const
virtual void Solve(int m, int n, real_t *X) const
virtual bool Factor(int m, real_t TOL=0.0)
Compute the Cholesky factorization of the current matrix.
void LSolve(int m, int n, real_t *X) const
void LMult(int m, int n, real_t *X) const
const Vector & Eigenvector(int i)
Definition densemat.hpp:907
DenseMatrixEigensystem(DenseMatrix &m)
DenseMatrix & Eigenvectors()
Definition densemat.hpp:905
DenseMatrixGeneralizedEigensystem(DenseMatrix &a, DenseMatrix &b, bool left_eigen_vectors=false, bool right_eigen_vectors=false)
void TestInversion()
Print the numerical conditioning of the inversion: ||A^{-1} A - I||.
int Size() const
Get the size of the inverse matrix.
Definition densemat.hpp:850
DenseMatrixInverse(bool spd_=false)
Default constructor.
Definition densemat.hpp:840
virtual ~DenseMatrixInverse()
Destroys dense inverse matrix.
void Factor()
Factor the current DenseMatrix, *a.
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication with the inverse of dense matrix.
real_t Det() const
Compute the determinant of the original DenseMatrix using the LU factors.
Definition densemat.hpp:878
void GetInverseMatrix(DenseMatrix &Ainv) const
Compute and return the inverse matrix in Ainv.
void Mult(DenseMatrix &X) const
Multiply the inverse matrix by another matrix: X <- A^{-1} X.
Definition densemat.hpp:870
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Class for Singular Value Decomposition of a DenseMatrix.
Definition densemat.hpp:956
Vector & Singularvalues()
Return singular values.
DenseMatrix & RightSingularvectors()
Return right singular vectors.
MFEM_DEPRECATED DenseMatrixSVD(DenseMatrix &M, bool left_singular_vectors=false, bool right_singular_vectors=false)
Constructor for the DenseMatrixSVD.
DenseMatrix & LeftSingularvectors()
Return left singular vectors.
real_t Singularvalue(int i)
Return specific singular value.
void Eval(DenseMatrix &M)
Evaluate the SVD.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void GetDiag(Vector &d) const
Returns the diagonal of the matrix.
void CopyMNDiag(real_t c, int n, int row_offset, int col_offset)
Copy c on the diagonal of size n to *this at row_offset, col_offset.
void AddMult_a(real_t a, const Vector &x, Vector &y) const
y += a * A.x
Definition densemat.cpp:346
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
Definition densemat.cpp:220
void AddMultTranspose_a(real_t a, const Vector &x, Vector &y) const
y += a * A^t x
Definition densemat.cpp:364
void Set(real_t alpha, const real_t *A)
Set the matrix to alpha * A, assuming that A has the same dimensions as the matrix and uses column-ma...
Definition densemat.cpp:619
void TestInversion()
Invert and print the numerical conditioning of the inversion.
void CopyExceptMN(const DenseMatrix &A, int m, int n)
Copy All rows and columns except m and n from A.
void CopyCols(const DenseMatrix &A, int col1, int col2)
Copy columns col1 through col2 from A to *this.
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
Definition densemat.cpp:262
void Transpose()
(*this) = (*this)^t
void AddToVector(int offset, Vector &v) const
Add the matrix 'data' to the Vector 'v' at the given 'offset'.
void Threshold(real_t eps)
Replace small entries, abs(a_ij) <= eps, with zero.
const real_t * HostRead() const
Shortcut for mfem::Read(GetMemory(), TotalSize(), false).
Definition densemat.hpp:470
void CalcEigenvalues(real_t *lambda, real_t *vec) const
void Eigenvalues(Vector &ev)
Compute eigenvalues of A x = ev x where A = *this.
Definition densemat.hpp:272
DenseMatrix(real_t *d, int h, int w)
Construct a DenseMatrix using an existing data array.
Definition densemat.hpp:58
real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev).
Definition densemat.hpp:482
void RightScaling(const Vector &s)
RightScaling: this = this * diag(s);.
Definition densemat.cpp:427
void Norm2(Vector &v) const
Take the 2-norm of the columns of A and store in v.
Definition densemat.hpp:256
real_t * HostWrite()
Shortcut for mfem::Write(GetMemory(), TotalSize(), false).
Definition densemat.hpp:478
virtual void PrintMatlab(std::ostream &out=mfem::out) const
Prints operator in Matlab format.
void SetRow(int r, const real_t *row)
void SymmetricScaling(const Vector &s)
SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
Definition densemat.cpp:456
void Getl1Diag(Vector &l) const
Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|.
void GetColumnReference(int c, Vector &col)
Definition densemat.hpp:312
real_t & operator()(int i, int j)
Returns reference to a_{ij}.
real_t InnerProduct(const real_t *x, const real_t *y) const
Compute y^t A x.
Definition densemat.cpp:383
void Eigenvalues(DenseMatrix &b, Vector &ev, DenseMatrix &evect)
Compute generalized eigenvalues of A x = ev B x, where A = *this.
Definition densemat.hpp:289
real_t * Data() const
Returns the matrix data array.
Definition densemat.hpp:111
const Memory< real_t > & GetMemory() const
Definition densemat.hpp:118
void GetSubMatrix(const Array< int > &idx, DenseMatrix &A) const
real_t * GetData() const
Returns the matrix data array.
Definition densemat.hpp:115
void Eigensystem(DenseMatrix &b, Vector &ev, DenseMatrix &evect)
Definition densemat.hpp:294
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
virtual ~DenseMatrix()
Destroys dense matrix.
virtual void AddMult(const Vector &x, Vector &y, const real_t a=1.0) const
y += a * A.x
Definition densemat.cpp:299
void Reset(real_t *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition densemat.hpp:90
virtual real_t & Elem(int i, int j)
Returns reference to a_{ij}.
Definition densemat.cpp:210
void SetCol(int c, const real_t *col)
void ClearExternalData()
Definition densemat.hpp:95
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
void Invert()
Replaces the current matrix with its inverse.
Definition densemat.cpp:726
real_t Weight() const
Definition densemat.cpp:592
DenseMatrix & operator+=(const real_t *m)
Definition densemat.cpp:681
void Neg()
(*this) = -(*this)
Definition densemat.cpp:717
real_t operator*(const DenseMatrix &m) const
Matrix inner product: tr(A^t B)
Definition densemat.cpp:247
void CopyMNt(const DenseMatrix &A, int row_offset, int col_offset)
Copy matrix A^t to the location in *this at row_offset, col_offset.
void AdjustDofDirection(Array< int > &dofs)
void Eigenvalues(Vector &ev, DenseMatrix &evect)
Compute eigenvalues and eigenvectors of A x = ev x where A = *this.
Definition densemat.hpp:276
void InvRightScaling(const Vector &s)
InvRightScaling: this = this * diag(1./s);.
Definition densemat.cpp:442
const real_t * GetColumn(int col) const
Definition densemat.hpp:310
void Eigensystem(Vector &ev, DenseMatrix &evect)
Compute eigenvalues and eigenvectors of A x = ev x where A = *this.
Definition densemat.hpp:280
real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
Definition densemat.hpp:474
virtual MatrixInverse * Inverse() const
Returns a pointer to the inverse matrix.
Definition densemat.cpp:530
real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false).
Definition densemat.hpp:486
void SingularValues(Vector &sv) const
real_t FNorm() const
Compute the Frobenius norm of the matrix.
Definition densemat.hpp:266
void InvLeftScaling(const Vector &s)
InvLeftScaling this = diag(1./s) * this.
Definition densemat.cpp:414
DenseMatrix(const T(&values)[M][N])
Definition densemat.hpp:64
virtual void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
void SetSubMatrix(const Array< int > &idx, const DenseMatrix &A)
Set (*this)(idx[i],idx[j]) = A(i,j)
virtual void PrintT(std::ostream &out=mfem::out, int width_=4) const
Prints the transpose matrix to stream out.
void AddSubMatrix(const Array< int > &idx, const DenseMatrix &A)
(*this)(idx[i],idx[j]) += A(i,j)
real_t Trace() const
Trace of a square matrix.
Definition densemat.cpp:511
real_t * GetColumn(int col)
Definition densemat.hpp:309
void Diag(real_t c, int n)
Creates n x n diagonal matrix with diagonal elements c.
void SquareRootInverse()
Replaces the current matrix with its square root inverse.
Definition densemat.cpp:853
virtual void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const
y += a * A^t x
Definition densemat.cpp:322
DenseMatrix & operator*=(real_t c)
Definition densemat.cpp:707
void Swap(DenseMatrix &other)
const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
Definition densemat.hpp:466
bool OwnsData() const
Return the DenseMatrix data (host pointer) ownership flag.
Definition densemat.hpp:121
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0<=i.
void GetRowSums(Vector &l) const
Compute the row sums of the DenseMatrix.
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
void InvSymmetricScaling(const Vector &s)
InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
Definition densemat.cpp:484
int Rank(real_t tol) const
void UseExternalData(real_t *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition densemat.hpp:80
void Add(const real_t c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition densemat.cpp:628
void Clear()
Delete the matrix data array (if owned) and reset the matrix state.
Definition densemat.hpp:98
void LeftScaling(const Vector &s)
LeftScaling this = diag(s) * this.
Definition densemat.cpp:401
DenseMatrix & operator-=(const DenseMatrix &m)
Definition densemat.cpp:694
real_t InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition densemat.hpp:198
std::size_t MemoryUsage() const
Definition densemat.hpp:463
Memory< real_t > & GetMemory()
Definition densemat.hpp:117
real_t CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
void GradToCurl(DenseMatrix &curl)
void GetColumn(int c, Vector &col) const
void GradToVectorCurl2D(DenseMatrix &curl)
DenseMatrix & operator=(real_t c)
Sets the matrix elements equal to constant c.
Definition densemat.cpp:648
void Eigenvalues(DenseMatrix &b, Vector &ev)
Definition densemat.hpp:285
real_t MaxMaxNorm() const
Compute the norm ||A|| = max_{ij} |A_{ij}|.
Definition densemat.cpp:911
void Norm2(real_t *v) const
Take the 2-norm of the columns of A and store in v.
Definition densemat.cpp:898
void GetFromVector(int offset, const Vector &v)
Get the matrix 'data' from the Vector 'v' at the given 'offset'.
int CheckFinite() const
Definition densemat.hpp:452
int Size() const
For backward compatibility define Size to be synonym of Width()
Definition densemat.hpp:102
real_t FNorm2() const
Compute the square of the Frobenius norm of the matrix.
Definition densemat.hpp:269
void GetRow(int r, Vector &row) const
void CopyRows(const DenseMatrix &A, int row1, int row2)
Copy rows row1 through row2 from A to *this.
real_t Det() const
Definition densemat.cpp:535
void GradToDiv(Vector &div)
void Set(real_t alpha, const DenseMatrix &A)
Set the matrix to alpha * A.
Definition densemat.hpp:220
Rank 3 tensor (array of matrices)
DenseTensor(int i, int j, int k, MemoryType mt)
real_t * HostWrite()
Shortcut for mfem::Write(GetMemory(), TotalSize(), false).
DenseTensor(const DenseTensor &other)
Copy constructor: deep copy.
void SetSize(int i, int j, int k, MemoryType mt_=MemoryType::PRESERVE)
real_t * GetData(int k)
real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false).
real_t & operator()(int i, int j, int k)
Memory< real_t > & GetMemory()
void UseExternalData(real_t *ext_data, int i, int j, int k)
int SizeJ() const
void AddMult(const Table &elem_dof, const Vector &x, Vector &y) const
int TotalSize() const
DenseTensor(real_t *d, int i, int j, int k)
const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
const real_t & operator()(int i, int j, int k) const
real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
std::size_t MemoryUsage() const
const real_t * Data() const
DenseTensor(int i, int j, int k)
real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev).
int SizeI() const
const real_t * HostRead() const
Shortcut for mfem::Read(GetMemory(), TotalSize(), false).
DenseTensor & operator=(real_t c)
Sets the tensor elements equal to constant c.
const DenseMatrix & operator()(int k) const
int SizeK() const
void Swap(DenseTensor &t)
const real_t * GetData(int k) const
DenseMatrix & operator()(int k)
const Memory< real_t > & GetMemory() const
virtual void GetInverseMatrix(int m, real_t *X) const
Definition densemat.hpp:652
virtual real_t Det(int m) const
Definition densemat.hpp:641
Factors(real_t *data_)
Definition densemat.hpp:633
virtual void Solve(int m, int n, real_t *X) const
Definition densemat.hpp:647
real_t * data
Definition densemat.hpp:629
virtual ~Factors()
Definition densemat.hpp:657
virtual bool Factor(int m, real_t TOL=0.0)
Definition densemat.hpp:635
A class to initialize the size of a Tensor.
Definition dtensor.hpp:55
void LSolve(int m, int n, real_t *X) const
static void SubMult(int m, int n, int r, const real_t *A21, const real_t *X1, real_t *X2)
virtual void GetInverseMatrix(int m, real_t *X) const
Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1}.
void Mult(int m, int n, real_t *X) const
void USolve(int m, int n, real_t *X) const
virtual bool Factor(int m, real_t TOL=0.0)
Compute the LU factorization of the current matrix.
void BlockFactor(int m, int n, real_t *A12, real_t *A21, real_t *A22) const
void BlockForwSolve(int m, int n, int r, const real_t *L21, real_t *B1, real_t *B2) const
void RightSolve(int m, int n, real_t *X) const
LUFactors(real_t *data_, int *ipiv_)
Definition densemat.hpp:677
virtual real_t Det(int m) const
virtual void Solve(int m, int n, real_t *X) const
void BlockBackSolve(int m, int n, int r, const real_t *U12, const real_t *X2, real_t *Y1) const
static const int ipiv_base
Definition densemat.hpp:668
Abstract data type for matrix inverse.
Definition matrix.hpp:63
Abstract data type matrix.
Definition matrix.hpp:28
Class used by MFEM to store pointers to host and/or device memory.
int Capacity() const
Return the size of the allocated memory.
bool OwnsHostPtr() const
Return true if the host pointer is owned. Ownership indicates whether the pointer will be deleted by ...
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
void CopyFrom(const Memory &src, int size)
Copy size entries from src to *this.
void Reset()
Reset the memory to be empty, ensuring that Delete() will be a no-op.
void Wrap(T *ptr, int size, bool own)
Wrap an externally allocated host pointer, ptr with the current host memory type returned by MemoryMa...
void Delete()
Delete the owned pointers and reset the Memory object.
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
Abstract operator.
Definition operator.hpp:25
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition operator.hpp:93
Vector data type.
Definition vector.hpp:80
void SetDataAndSize(real_t *d, int s)
Set the Vector data and size.
Definition vector.hpp:175
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
void SetData(real_t *d)
Definition vector.hpp:168
Vector beta
const real_t alpha
Definition ex15.cpp:369
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
void CalcOrtho(const DenseMatrix &J, Vector &n)
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true,...
Definition device.hpp:320
void Swap(Array< T > &, Array< T > &)
Definition array.hpp:648
void AddMult_a_ABt(real_t a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
void mfem_error(const char *msg)
Definition error.cpp:154
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:476
void BatchLUSolve(const DenseTensor &Mlu, const Array< int > &P, Vector &X)
Solve batch linear systems.
void AddMult_a(real_t alpha, const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += alpha * B * C.
void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat)
Calculate the transposed adjugate of a matrix (for NxN matrices, N=1,2,3)
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
void MultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
void AddMultVVt(const Vector &v, DenseMatrix &VVt)
VVt += v v^t.
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true,...
Definition device.hpp:337
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt = A D A^t, where D is diagonal.
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
T * ReadWrite(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read+write access to mem with the mfem::Device's DeviceMemoryClass,...
Definition device.hpp:354
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
void AddMult_a_VWt(const real_t a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
void AddMult_a_VVt(const real_t a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
void BatchLUFactor(DenseTensor &Mlu, Array< int > &P, const real_t TOL)
Compute the LU factorization of a batch of matrices.
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
int CheckFinite(const real_t *v, const int n)
Definition vector.hpp:518
void AddMult_a_AtB(real_t a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
AtB += a * A^t * B.
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
bool LinearSolve(DenseMatrix &A, real_t *X, real_t TOL)
Solves the dense linear system, A * X = B for X
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
void AddMult_a_AAt(real_t a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva)
Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
void AddMultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
float real_t
Definition config.hpp:43
void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += B * C.
void AddMultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
AtB += A^t * B.
MemoryType
Memory types supported by MFEM.
void Mult_a_AAt(real_t a, const DenseMatrix &A, DenseMatrix &AAt)
AAt = a * A * A^t.
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
RefCoord t[3]
RefCoord s[3]