MFEM
v4.1.0
Finite element discretization library
|
Data type dense matrix using column-major storage. More...
#include <densemat.hpp>
Public Member Functions | |
DenseMatrix () | |
DenseMatrix (const DenseMatrix &) | |
Copy constructor. More... | |
DenseMatrix (int s) | |
Creates square matrix of size s. More... | |
DenseMatrix (int m, int n) | |
Creates rectangular matrix of size m x n. More... | |
DenseMatrix (const DenseMatrix &mat, char ch) | |
Creates rectangular matrix equal to the transpose of mat. More... | |
DenseMatrix (double *d, int h, int w) | |
Construct a DenseMatrix using an existing data array. More... | |
void | UseExternalData (double *d, int h, int w) |
Change the data array and the size of the DenseMatrix. More... | |
void | Reset (double *d, int h, int w) |
Change the data array and the size of the DenseMatrix. More... | |
void | ClearExternalData () |
void | Clear () |
Delete the matrix data array (if owned) and reset the matrix state. More... | |
int | Size () const |
For backward compatibility define Size to be synonym of Width() More... | |
void | SetSize (int s) |
Change the size of the DenseMatrix to s x s. More... | |
void | SetSize (int h, int w) |
Change the size of the DenseMatrix to h x w. More... | |
double * | Data () const |
Returns the matrix data array. More... | |
double * | GetData () const |
Returns the matrix data array. More... | |
Memory< double > & | GetMemory () |
const Memory< double > & | GetMemory () const |
bool | OwnsData () const |
Return the DenseMatrix data (host pointer) ownership flag. More... | |
double & | operator() (int i, int j) |
Returns reference to a_{ij}. More... | |
const double & | operator() (int i, int j) const |
Returns constant reference to a_{ij}. More... | |
double | operator* (const DenseMatrix &m) const |
Matrix inner product: tr(A^t B) More... | |
double | Trace () const |
Trace of a square matrix. More... | |
virtual double & | Elem (int i, int j) |
Returns reference to a_{ij}. More... | |
virtual const double & | Elem (int i, int j) const |
Returns constant reference to a_{ij}. More... | |
void | Mult (const double *x, double *y) const |
Matrix vector multiplication. More... | |
virtual void | Mult (const Vector &x, Vector &y) const |
Matrix vector multiplication. More... | |
void | MultTranspose (const double *x, double *y) const |
Multiply a vector with the transpose matrix. More... | |
virtual void | MultTranspose (const Vector &x, Vector &y) const |
Multiply a vector with the transpose matrix. More... | |
void | AddMult (const Vector &x, Vector &y) const |
y += A.x More... | |
void | AddMultTranspose (const Vector &x, Vector &y) const |
y += A^t x More... | |
void | AddMult_a (double a, const Vector &x, Vector &y) const |
y += a * A.x More... | |
void | AddMultTranspose_a (double a, const Vector &x, Vector &y) const |
y += a * A^t x More... | |
double | InnerProduct (const double *x, const double *y) const |
Compute y^t A x. More... | |
void | LeftScaling (const Vector &s) |
LeftScaling this = diag(s) * this. More... | |
void | InvLeftScaling (const Vector &s) |
InvLeftScaling this = diag(1./s) * this. More... | |
void | RightScaling (const Vector &s) |
RightScaling: this = this * diag(s);. More... | |
void | InvRightScaling (const Vector &s) |
InvRightScaling: this = this * diag(1./s);. More... | |
void | SymmetricScaling (const Vector &s) |
SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s)) More... | |
void | InvSymmetricScaling (const Vector &s) |
InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s)) More... | |
double | InnerProduct (const Vector &x, const Vector &y) const |
Compute y^t A x. More... | |
virtual MatrixInverse * | Inverse () const |
Returns a pointer to the inverse matrix. More... | |
void | Invert () |
Replaces the current matrix with its inverse. More... | |
void | SquareRootInverse () |
Replaces the current matrix with its square root inverse. More... | |
double | Det () const |
double | Weight () const |
void | Set (double alpha, const double *A) |
Set the matrix to alpha * A, assuming that A has the same dimensions as the matrix and uses column-major layout. More... | |
void | Set (double alpha, const DenseMatrix &A) |
Set the matrix to alpha * A. More... | |
void | Add (const double c, const DenseMatrix &A) |
Adds the matrix A multiplied by the number c to the matrix. More... | |
DenseMatrix & | operator= (double c) |
Sets the matrix elements equal to constant c. More... | |
DenseMatrix & | operator= (const double *d) |
Copy the matrix entries from the given array. More... | |
DenseMatrix & | operator= (const DenseMatrix &m) |
Sets the matrix size and elements equal to those of m. More... | |
DenseMatrix & | operator+= (const double *m) |
DenseMatrix & | operator+= (const DenseMatrix &m) |
DenseMatrix & | operator-= (const DenseMatrix &m) |
DenseMatrix & | operator*= (double c) |
void | Neg () |
(*this) = -(*this) More... | |
void | Norm2 (double *v) const |
Take the 2-norm of the columns of A and store in v. More... | |
double | MaxMaxNorm () const |
Compute the norm ||A|| = max_{ij} |A_{ij}|. More... | |
double | FNorm () const |
Compute the Frobenius norm of the matrix. More... | |
double | FNorm2 () const |
Compute the square of the Frobenius norm of the matrix. More... | |
void | Eigenvalues (Vector &ev) |
Compute eigenvalues of A x = ev x where A = *this. More... | |
void | Eigenvalues (Vector &ev, DenseMatrix &evect) |
Compute eigenvalues and eigenvectors of A x = ev x where A = *this. More... | |
void | Eigensystem (Vector &ev, DenseMatrix &evect) |
Compute eigenvalues and eigenvectors of A x = ev x where A = *this. More... | |
void | Eigenvalues (DenseMatrix &b, Vector &ev) |
void | Eigenvalues (DenseMatrix &b, Vector &ev, DenseMatrix &evect) |
Compute generalized eigenvalues of A x = ev B x, where A = *this. More... | |
void | Eigensystem (DenseMatrix &b, Vector &ev, DenseMatrix &evect) |
void | SingularValues (Vector &sv) const |
int | Rank (double tol) const |
double | CalcSingularvalue (const int i) const |
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3. More... | |
void | CalcEigenvalues (double *lambda, double *vec) const |
void | GetRow (int r, Vector &row) const |
void | GetColumn (int c, Vector &col) const |
double * | GetColumn (int col) |
const double * | GetColumn (int col) const |
void | GetColumnReference (int c, Vector &col) |
void | SetRow (int r, const double *row) |
void | SetRow (int r, const Vector &row) |
void | SetCol (int c, const double *col) |
void | SetCol (int c, const Vector &col) |
void | SetRow (int row, double value) |
Set all entries of a row to the specified value. More... | |
void | SetCol (int col, double value) |
Set all entries of a column to the specified value. More... | |
void | GetDiag (Vector &d) const |
Returns the diagonal of the matrix. More... | |
void | Getl1Diag (Vector &l) const |
Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|. More... | |
void | GetRowSums (Vector &l) const |
Compute the row sums of the DenseMatrix. More... | |
void | Diag (double c, int n) |
Creates n x n diagonal matrix with diagonal elements c. More... | |
void | Diag (double *diag, int n) |
Creates n x n diagonal matrix with diagonal given by diag. More... | |
void | Transpose () |
(*this) = (*this)^t More... | |
void | Transpose (const DenseMatrix &A) |
(*this) = A^t More... | |
void | Symmetrize () |
(*this) = 1/2 ((*this) + (*this)^t) More... | |
void | Lump () |
void | GradToCurl (DenseMatrix &curl) |
void | GradToDiv (Vector &div) |
void | CopyRows (const DenseMatrix &A, int row1, int row2) |
Copy rows row1 through row2 from A to *this. More... | |
void | CopyCols (const DenseMatrix &A, int col1, int col2) |
Copy columns col1 through col2 from A to *this. More... | |
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. More... | |
void | CopyMN (const DenseMatrix &A, int row_offset, int col_offset) |
Copy matrix A to the location in *this at row_offset, col_offset. More... | |
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. More... | |
void | CopyMN (const DenseMatrix &A, int m, int n, int Aro, int Aco, int row_offset, int col_offset) |
void | CopyMNDiag (double c, int n, int row_offset, int col_offset) |
Copy c on the diagonal of size n to *this at row_offset, col_offset. More... | |
void | CopyMNDiag (double *diag, int n, int row_offset, int col_offset) |
Copy diag on the diagonal of size n to *this at row_offset, col_offset. More... | |
void | CopyExceptMN (const DenseMatrix &A, int m, int n) |
Copy All rows and columns except m and n from A. More... | |
void | AddMatrix (DenseMatrix &A, int ro, int co) |
Perform (ro+i,co+j)+=A(i,j) for 0<=i<A.Height, 0<=j<A.Width. More... | |
void | AddMatrix (double a, const DenseMatrix &A, int ro, int co) |
Perform (ro+i,co+j)+=a*A(i,j) for 0<=i<A.Height, 0<=j<A.Width. More... | |
void | AddToVector (int offset, Vector &v) const |
Add the matrix 'data' to the Vector 'v' at the given 'offset'. More... | |
void | GetFromVector (int offset, const Vector &v) |
Get the matrix 'data' from the Vector 'v' at the given 'offset'. More... | |
void | AdjustDofDirection (Array< int > &dofs) |
void | Threshold (double eps) |
Replace small entries, abs(a_ij) <= eps, with zero. More... | |
int | CheckFinite () const |
virtual void | Print (std::ostream &out=mfem::out, int width_=4) const |
Prints matrix to stream out. More... | |
virtual void | PrintMatlab (std::ostream &out=mfem::out) const |
virtual void | PrintT (std::ostream &out=mfem::out, int width_=4) const |
Prints the transpose matrix to stream out. More... | |
void | TestInversion () |
Invert and print the numerical conditioning of the inversion. More... | |
long | MemoryUsage () const |
const double * | Read (bool on_dev=true) const |
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev). More... | |
const double * | HostRead () const |
Shortcut for mfem::Read(GetMemory(), TotalSize(), false). More... | |
double * | Write (bool on_dev=true) |
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev). More... | |
double * | HostWrite () |
Shortcut for mfem::Write(GetMemory(), TotalSize(), false). More... | |
double * | ReadWrite (bool on_dev=true) |
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev). More... | |
double * | HostReadWrite () |
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false). More... | |
virtual | ~DenseMatrix () |
Destroys dense matrix. More... | |
Public Member Functions inherited from mfem::Matrix | |
Matrix (int s) | |
Creates a square matrix of size s. More... | |
Matrix (int h, int w) | |
Creates a matrix of the given height and width. More... | |
bool | IsSquare () const |
Returns whether the matrix is a square matrix. More... | |
virtual void | Finalize (int) |
Finalizes the matrix initialization. More... | |
virtual | ~Matrix () |
Destroys matrix. More... | |
Public Member Functions inherited from mfem::Operator | |
void | InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const |
Initializes memory for true vectors of linear system. More... | |
Operator (int s=0) | |
Construct a square Operator with given size s (default 0). More... | |
Operator (int h, int w) | |
Construct an Operator with the given height (output size) and width (input size). More... | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows(). More... | |
int | NumRows () const |
Get the number of rows (size of output) of the Operator. Synonym with Height(). More... | |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols(). More... | |
int | NumCols () const |
Get the number of columns (size of input) of the Operator. Synonym with Width(). More... | |
virtual MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. More... | |
virtual Operator & | GetGradient (const Vector &x) const |
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More... | |
virtual const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
virtual const Operator * | GetOutputProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetOutputRestriction () const |
Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0) |
Form a constrained linear system using a matrix-free approach. More... | |
void | FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B) |
Form a column-constrained linear system using a matrix-free approach. More... | |
virtual void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem() or Operator::FormRectangularLinearSystem(). More... | |
void | FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this square operator. More... | |
void | FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator (including constraints). More... | |
void | FormDiscreteOperator (Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator. More... | |
void | PrintMatlab (std::ostream &out, int n=0, int m=0) const |
Prints operator with input size n and output size m in Matlab format. More... | |
virtual | ~Operator () |
Virtual destructor. More... | |
Type | GetType () const |
Return the type ID of the Operator class. More... | |
Friends | |
class | DenseTensor |
class | DenseMatrixInverse |
Additional Inherited Members | |
Public Types inherited from mfem::Matrix | |
enum | DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP } |
Public Types inherited from mfem::Operator | |
enum | Type { ANY_TYPE, MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ, PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE, PETSC_MATGENERIC, Complex_Operator, MFEM_ComplexSparseMat, Complex_Hypre_ParCSR } |
Enumeration defining IDs for some classes derived from Operator. More... | |
Protected Member Functions inherited from mfem::Operator | |
void | FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout) |
see FormSystemOperator() More... | |
void | FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout) |
see FormRectangularSystemOperator() More... | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, taking in input/output Prolongation matrices. More... | |
Protected Attributes inherited from mfem::Operator | |
int | height |
Dimension of the output / number of rows in the matrix. More... | |
int | width |
Dimension of the input / number of columns in the matrix. More... | |
Data type dense matrix using column-major storage.
Definition at line 23 of file densemat.hpp.
mfem::DenseMatrix::DenseMatrix | ( | ) |
Default constructor for DenseMatrix. Sets data = NULL and height = width = 0.
Definition at line 71 of file densemat.cpp.
mfem::DenseMatrix::DenseMatrix | ( | const DenseMatrix & | m | ) |
Copy constructor.
Definition at line 76 of file densemat.cpp.
|
explicit |
Creates square matrix of size s.
Definition at line 91 of file densemat.cpp.
mfem::DenseMatrix::DenseMatrix | ( | int | m, |
int | n | ||
) |
Creates rectangular matrix of size m x n.
Definition at line 105 of file densemat.cpp.
mfem::DenseMatrix::DenseMatrix | ( | const DenseMatrix & | mat, |
char | ch | ||
) |
Creates rectangular matrix equal to the transpose of mat.
Definition at line 121 of file densemat.cpp.
|
inline |
Construct a DenseMatrix using an existing data array.
The DenseMatrix does not assume ownership of the data array, i.e. it will not delete the array.
Definition at line 58 of file densemat.hpp.
|
virtual |
Destroys dense matrix.
Definition at line 1919 of file densemat.cpp.
void mfem::DenseMatrix::Add | ( | const double | c, |
const DenseMatrix & | A | ||
) |
Adds the matrix A multiplied by the number c to the matrix.
Definition at line 542 of file densemat.cpp.
void mfem::DenseMatrix::AddMatrix | ( | DenseMatrix & | A, |
int | ro, | ||
int | co | ||
) |
Perform (ro+i,co+j)+=A(i,j) for 0<=i<A.Height, 0<=j<A.Width.
Definition at line 1665 of file densemat.cpp.
void mfem::DenseMatrix::AddMatrix | ( | double | a, |
const DenseMatrix & | A, | ||
int | ro, | ||
int | co | ||
) |
Perform (ro+i,co+j)+=a*A(i,j) for 0<=i<A.Height, 0<=j<A.Width.
Definition at line 1695 of file densemat.cpp.
y += A.x
Definition at line 224 of file densemat.cpp.
y += a * A.x
Definition at line 260 of file densemat.cpp.
y += A^t x
Definition at line 242 of file densemat.cpp.
y += a * A^t x
Definition at line 278 of file densemat.cpp.
void mfem::DenseMatrix::AddToVector | ( | int | offset, |
Vector & | v | ||
) | const |
Add the matrix 'data' to the Vector 'v' at the given 'offset'.
Definition at line 1725 of file densemat.cpp.
void mfem::DenseMatrix::AdjustDofDirection | ( | Array< int > & | dofs | ) |
If (dofs[i] < 0 and dofs[j] >= 0) or (dofs[i] >= 0 and dofs[j] < 0) then (*this)(i,j) = -(*this)(i,j).
Definition at line 1747 of file densemat.cpp.
void mfem::DenseMatrix::CalcEigenvalues | ( | double * | lambda, |
double * | vec | ||
) | const |
Return the eigenvalues (in increasing order) and eigenvectors of a 2x2 or 3x3 symmetric matrix.
Definition at line 1249 of file densemat.cpp.
double mfem::DenseMatrix::CalcSingularvalue | ( | const int | i | ) | const |
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition at line 1224 of file densemat.cpp.
|
inline |
Count the number of entries in the matrix for which isfinite is false, i.e. the entry is a NaN or +/-Inf.
Definition at line 356 of file densemat.hpp.
|
inline |
Delete the matrix data array (if owned) and reset the matrix state.
Definition at line 83 of file densemat.hpp.
|
inline |
Clear the data array and the dimensions of the DenseMatrix. This method should not be used with DenseMatrix that owns its current data array.
Definition at line 80 of file densemat.hpp.
void mfem::DenseMatrix::CopyCols | ( | const DenseMatrix & | A, |
int | col1, | ||
int | col2 | ||
) |
Copy columns col1 through col2 from A to *this.
Definition at line 1521 of file densemat.cpp.
void mfem::DenseMatrix::CopyExceptMN | ( | const DenseMatrix & | A, |
int | m, | ||
int | n | ||
) |
Copy All rows and columns except m and n from A.
Definition at line 1639 of file densemat.cpp.
void mfem::DenseMatrix::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.
Definition at line 1534 of file densemat.cpp.
void mfem::DenseMatrix::CopyMN | ( | const DenseMatrix & | A, |
int | row_offset, | ||
int | col_offset | ||
) |
Copy matrix A to the location in *this at row_offset, col_offset.
Definition at line 1547 of file densemat.cpp.
void mfem::DenseMatrix::CopyMN | ( | const DenseMatrix & | A, |
int | m, | ||
int | n, | ||
int | Aro, | ||
int | Aco, | ||
int | row_offset, | ||
int | col_offset | ||
) |
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this at row_offset, col_offset
Definition at line 1573 of file densemat.cpp.
void mfem::DenseMatrix::CopyMNDiag | ( | double | c, |
int | n, | ||
int | row_offset, | ||
int | col_offset | ||
) |
Copy c on the diagonal of size n to *this at row_offset, col_offset.
Definition at line 1604 of file densemat.cpp.
void mfem::DenseMatrix::CopyMNDiag | ( | double * | diag, |
int | n, | ||
int | row_offset, | ||
int | col_offset | ||
) |
Copy diag on the diagonal of size n to *this at row_offset, col_offset.
Definition at line 1621 of file densemat.cpp.
void mfem::DenseMatrix::CopyMNt | ( | const DenseMatrix & | A, |
int | row_offset, | ||
int | col_offset | ||
) |
Copy matrix A^t to the location in *this at row_offset, col_offset.
Definition at line 1560 of file densemat.cpp.
void mfem::DenseMatrix::CopyRows | ( | const DenseMatrix & | A, |
int | row1, | ||
int | row2 | ||
) |
Copy rows row1 through row2 from A to *this.
Definition at line 1508 of file densemat.cpp.
|
inline |
Returns the matrix data array.
Definition at line 96 of file densemat.hpp.
double mfem::DenseMatrix::Det | ( | ) | const |
Calculates the determinant of the matrix (optimized for 2x2, 3x3, and 4x4 matrices)
Definition at line 449 of file densemat.cpp.
void mfem::DenseMatrix::Diag | ( | double | c, |
int | n | ||
) |
Creates n x n diagonal matrix with diagonal elements c.
Definition at line 1346 of file densemat.cpp.
void mfem::DenseMatrix::Diag | ( | double * | diag, |
int | n | ||
) |
Creates n x n diagonal matrix with diagonal given by diag.
Definition at line 1361 of file densemat.cpp.
|
inline |
Compute eigenvalues and eigenvectors of A x = ev x where A = *this.
Definition at line 238 of file densemat.hpp.
|
inline |
Compute generalized eigenvalues and eigenvectors of A x = ev B x, where A = *this
Definition at line 252 of file densemat.hpp.
|
inline |
Compute eigenvalues of A x = ev x where A = *this.
Definition at line 230 of file densemat.hpp.
|
inline |
Compute eigenvalues and eigenvectors of A x = ev x where A = *this.
Definition at line 234 of file densemat.hpp.
|
inline |
Compute generalized eigenvalues and eigenvectors of A x = ev B x, where A = *this
Definition at line 243 of file densemat.hpp.
|
inline |
Compute generalized eigenvalues of A x = ev B x, where A = *this.
Definition at line 247 of file densemat.hpp.
|
virtual |
|
virtual |
Returns constant reference to a_{ij}.
Implements mfem::Matrix.
Definition at line 168 of file densemat.cpp.
|
inline |
Compute the Frobenius norm of the matrix.
Definition at line 224 of file densemat.hpp.
|
inline |
Compute the square of the Frobenius norm of the matrix.
Definition at line 227 of file densemat.hpp.
void mfem::DenseMatrix::GetColumn | ( | int | c, |
Vector & | col | ||
) | const |
Definition at line 1287 of file densemat.cpp.
|
inline |
Definition at line 267 of file densemat.hpp.
|
inline |
Definition at line 268 of file densemat.hpp.
|
inline |
Definition at line 270 of file densemat.hpp.
|
inline |
Returns the matrix data array.
Definition at line 100 of file densemat.hpp.
void mfem::DenseMatrix::GetDiag | ( | Vector & | d | ) | const |
Returns the diagonal of the matrix.
Definition at line 1301 of file densemat.cpp.
void mfem::DenseMatrix::GetFromVector | ( | int | offset, |
const Vector & | v | ||
) |
Get the matrix 'data' from the Vector 'v' at the given 'offset'.
Definition at line 1736 of file densemat.cpp.
void mfem::DenseMatrix::Getl1Diag | ( | Vector & | l | ) | const |
Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|.
Definition at line 1315 of file densemat.cpp.
|
inline |
Definition at line 102 of file densemat.hpp.
|
inline |
Definition at line 103 of file densemat.hpp.
void mfem::DenseMatrix::GetRow | ( | int | r, |
Vector & | row | ||
) | const |
Definition at line 1271 of file densemat.cpp.
void mfem::DenseMatrix::GetRowSums | ( | Vector & | l | ) | const |
Compute the row sums of the DenseMatrix.
Definition at line 1332 of file densemat.cpp.
void mfem::DenseMatrix::GradToCurl | ( | DenseMatrix & | curl | ) |
Given a DShape matrix (from a scalar FE), stored in this, returns the CurlShape matrix. If *this is a N by D matrix, then curl is a D*N by D(D-1)/2 matrix. The size of curl must be set outside. The dimension D can be either 2 or 3.
Definition at line 1434 of file densemat.cpp.
void mfem::DenseMatrix::GradToDiv | ( | Vector & | div | ) |
Given a DShape matrix (from a scalar FE), stored in *this, returns the DivShape vector. If *this is a N by dim matrix, then div is a dim*N vector. The size of div must be set outside.
Definition at line 1493 of file densemat.cpp.
|
inline |
Shortcut for mfem::Read(GetMemory(), TotalSize(), false).
Definition at line 374 of file densemat.hpp.
|
inline |
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false).
Definition at line 390 of file densemat.hpp.
|
inline |
Shortcut for mfem::Write(GetMemory(), TotalSize(), false).
Definition at line 382 of file densemat.hpp.
double mfem::DenseMatrix::InnerProduct | ( | const double * | x, |
const double * | y | ||
) | const |
Compute y^t A x.
Definition at line 297 of file densemat.cpp.
Compute y^t A x.
Definition at line 167 of file densemat.hpp.
|
virtual |
Returns a pointer to the inverse matrix.
Implements mfem::Matrix.
Definition at line 444 of file densemat.cpp.
void mfem::DenseMatrix::Invert | ( | ) |
Replaces the current matrix with its inverse.
Definition at line 635 of file densemat.cpp.
void mfem::DenseMatrix::InvLeftScaling | ( | const Vector & | s | ) |
InvLeftScaling this = diag(1./s) * this.
Definition at line 328 of file densemat.cpp.
void mfem::DenseMatrix::InvRightScaling | ( | const Vector & | s | ) |
InvRightScaling: this = this * diag(1./s);.
Definition at line 356 of file densemat.cpp.
void mfem::DenseMatrix::InvSymmetricScaling | ( | const Vector & | s | ) |
InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
Definition at line 398 of file densemat.cpp.
void mfem::DenseMatrix::LeftScaling | ( | const Vector & | s | ) |
LeftScaling this = diag(s) * this.
Definition at line 315 of file densemat.cpp.
void mfem::DenseMatrix::Lump | ( | ) |
Definition at line 1420 of file densemat.cpp.
double mfem::DenseMatrix::MaxMaxNorm | ( | ) | const |
Compute the norm ||A|| = max_{ij} |A_{ij}|.
Definition at line 803 of file densemat.cpp.
|
inline |
Definition at line 367 of file densemat.hpp.
void mfem::DenseMatrix::Mult | ( | const double * | x, |
double * | y | ||
) | const |
Matrix vector multiplication.
Definition at line 173 of file densemat.cpp.
Matrix vector multiplication.
Implements mfem::Operator.
Definition at line 178 of file densemat.cpp.
void mfem::DenseMatrix::MultTranspose | ( | const double * | x, |
double * | y | ||
) | const |
Multiply a vector with the transpose matrix.
Definition at line 201 of file densemat.cpp.
Multiply a vector with the transpose matrix.
Reimplemented from mfem::Operator.
Definition at line 216 of file densemat.cpp.
void mfem::DenseMatrix::Neg | ( | ) |
(*this) = -(*this)
Definition at line 626 of file densemat.cpp.
void mfem::DenseMatrix::Norm2 | ( | double * | v | ) | const |
Take the 2-norm of the columns of A and store in v.
Definition at line 790 of file densemat.cpp.
|
inline |
Returns reference to a_{ij}.
Definition at line 865 of file densemat.hpp.
|
inline |
Returns constant reference to a_{ij}.
Definition at line 871 of file densemat.hpp.
double mfem::DenseMatrix::operator* | ( | const DenseMatrix & | m | ) | const |
Matrix inner product: tr(A^t B)
Definition at line 186 of file densemat.cpp.
DenseMatrix & mfem::DenseMatrix::operator*= | ( | double | c | ) |
Definition at line 616 of file densemat.cpp.
DenseMatrix & mfem::DenseMatrix::operator+= | ( | const double * | m | ) |
Definition at line 586 of file densemat.cpp.
DenseMatrix & mfem::DenseMatrix::operator+= | ( | const DenseMatrix & | m | ) |
Definition at line 596 of file densemat.cpp.
DenseMatrix & mfem::DenseMatrix::operator-= | ( | const DenseMatrix & | m | ) |
Definition at line 603 of file densemat.cpp.
DenseMatrix & mfem::DenseMatrix::operator= | ( | double | c | ) |
Sets the matrix elements equal to constant c.
Definition at line 553 of file densemat.cpp.
DenseMatrix & mfem::DenseMatrix::operator= | ( | const double * | d | ) |
Copy the matrix entries from the given array.
Definition at line 563 of file densemat.cpp.
DenseMatrix & mfem::DenseMatrix::operator= | ( | const DenseMatrix & | m | ) |
Sets the matrix size and elements equal to those of m.
Definition at line 573 of file densemat.cpp.
|
inline |
Return the DenseMatrix data (host pointer) ownership flag.
Definition at line 106 of file densemat.hpp.
|
virtual |
Prints matrix to stream out.
Reimplemented from mfem::Matrix.
Definition at line 1834 of file densemat.cpp.
|
virtual |
Definition at line 1860 of file densemat.cpp.
|
virtual |
Prints the transpose matrix to stream out.
Definition at line 1879 of file densemat.cpp.
int mfem::DenseMatrix::Rank | ( | double | tol | ) | const |
Definition at line 1209 of file densemat.cpp.
|
inline |
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
Definition at line 370 of file densemat.hpp.
|
inline |
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev).
Definition at line 386 of file densemat.hpp.
|
inline |
Change the data array and the size of the DenseMatrix.
The DenseMatrix does not assume ownership of the data array, i.e. it will not delete the new array d. This method will delete the current data array, if owned.
Definition at line 75 of file densemat.hpp.
void mfem::DenseMatrix::RightScaling | ( | const Vector & | s | ) |
RightScaling: this = this * diag(s);.
Definition at line 341 of file densemat.cpp.
void mfem::DenseMatrix::Set | ( | double | alpha, |
const double * | A | ||
) |
Set the matrix to alpha * A, assuming that A has the same dimensions as the matrix and uses column-major layout.
Definition at line 533 of file densemat.cpp.
|
inline |
Set the matrix to alpha * A.
Definition at line 189 of file densemat.hpp.
void mfem::DenseMatrix::SetCol | ( | int | c, |
const double * | col | ||
) |
Definition at line 1805 of file densemat.cpp.
void mfem::DenseMatrix::SetCol | ( | int | c, |
const Vector & | col | ||
) |
Definition at line 1814 of file densemat.cpp.
void mfem::DenseMatrix::SetCol | ( | int | col, |
double | value | ||
) |
Set all entries of a column to the specified value.
Definition at line 1782 of file densemat.cpp.
void mfem::DenseMatrix::SetRow | ( | int | r, |
const double * | row | ||
) |
Definition at line 1790 of file densemat.cpp.
void mfem::DenseMatrix::SetRow | ( | int | r, |
const Vector & | row | ||
) |
Definition at line 1799 of file densemat.cpp.
void mfem::DenseMatrix::SetRow | ( | int | row, |
double | value | ||
) |
Set all entries of a row to the specified value.
Definition at line 1774 of file densemat.cpp.
|
inline |
Change the size of the DenseMatrix to s x s.
Definition at line 90 of file densemat.hpp.
void mfem::DenseMatrix::SetSize | ( | int | h, |
int | w | ||
) |
Change the size of the DenseMatrix to h x w.
Definition at line 144 of file densemat.cpp.
void mfem::DenseMatrix::SingularValues | ( | Vector & | sv | ) | const |
Definition at line 1169 of file densemat.cpp.
|
inline |
For backward compatibility define Size to be synonym of Width()
Definition at line 87 of file densemat.hpp.
void mfem::DenseMatrix::SquareRootInverse | ( | ) |
Replaces the current matrix with its square root inverse.
Definition at line 745 of file densemat.cpp.
void mfem::DenseMatrix::SymmetricScaling | ( | const Vector & | s | ) |
SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
Definition at line 370 of file densemat.cpp.
void mfem::DenseMatrix::Symmetrize | ( | ) |
(*this) = 1/2 ((*this) + (*this)^t)
Definition at line 1409 of file densemat.cpp.
void mfem::DenseMatrix::TestInversion | ( | ) |
Invert and print the numerical conditioning of the inversion.
Definition at line 1905 of file densemat.cpp.
void mfem::DenseMatrix::Threshold | ( | double | eps | ) |
Replace small entries, abs(a_ij) <= eps, with zero.
Definition at line 1820 of file densemat.cpp.
double mfem::DenseMatrix::Trace | ( | ) | const |
Trace of a square matrix.
Definition at line 425 of file densemat.cpp.
void mfem::DenseMatrix::Transpose | ( | ) |
(*this) = (*this)^t
Definition at line 1376 of file densemat.cpp.
void mfem::DenseMatrix::Transpose | ( | const DenseMatrix & | A | ) |
(*this) = A^t
Definition at line 1398 of file densemat.cpp.
|
inline |
Change the data array and the size of the DenseMatrix.
The DenseMatrix does not assume ownership of the data array, i.e. it will not delete the data array d. This method should not be used with DenseMatrix that owns its current data array.
Definition at line 65 of file densemat.hpp.
double mfem::DenseMatrix::Weight | ( | ) | const |
Definition at line 506 of file densemat.cpp.
|
inline |
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
Definition at line 378 of file densemat.hpp.
|
friend |
Definition at line 26 of file densemat.hpp.
|
friend |
Definition at line 25 of file densemat.hpp.