MFEM  v4.6.0
Finite element discretization library
Public Member Functions | Friends | List of all members
mfem::DenseMatrix Class Reference

Data type dense matrix using column-major storage. More...

#include <densemat.hpp>

Inheritance diagram for mfem::DenseMatrix:
[legend]
Collaboration diagram for mfem::DenseMatrix:
[legend]

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...
 
template<int M, int N>
 DenseMatrix (const double(&values)[M][N])
 
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...
 
void Mult (const double *x, Vector &y) const
 Matrix vector multiplication. More...
 
void Mult (const Vector &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...
 
void MultTranspose (const double *x, Vector &y) const
 Multiply a vector with the transpose matrix. More...
 
void MultTranspose (const Vector &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...
 
virtual void AddMult (const Vector &x, Vector &y, const double a=1.0) const
 y += a * A.x More...
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const
 y += a * 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 MatrixInverseInverse () 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...
 
void Add (const double c, const double *A)
 
DenseMatrixoperator= (double c)
 Sets the matrix elements equal to constant c. More...
 
DenseMatrixoperator= (const double *d)
 Copy the matrix entries from the given array. More...
 
DenseMatrixoperator= (const DenseMatrix &m)
 Sets the matrix size and elements equal to those of m. More...
 
DenseMatrixoperator+= (const double *m)
 
DenseMatrixoperator+= (const DenseMatrix &m)
 
DenseMatrixoperator-= (const DenseMatrix &m)
 
DenseMatrixoperator*= (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...
 
void Norm2 (Vector &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 GradToVectorCurl2D (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 GetSubMatrix (const Array< int > &idx, DenseMatrix &A) const
 
void GetSubMatrix (const Array< int > &idx_i, const Array< int > &idx_j, DenseMatrix &A) const
 
void GetSubMatrix (int ibeg, int iend, DenseMatrix &A)
 
void GetSubMatrix (int ibeg, int iend, int jbeg, int jend, DenseMatrix &A)
 
void SetSubMatrix (const Array< int > &idx, const DenseMatrix &A)
 Set (*this)(idx[i],idx[j]) = A(i,j) More...
 
void SetSubMatrix (const Array< int > &idx_i, const Array< int > &idx_j, const DenseMatrix &A)
 Set (*this)(idx_i[i],idx_j[j]) = A(i,j) More...
 
void SetSubMatrix (int ibeg, const DenseMatrix &A)
 
void SetSubMatrix (int ibeg, int jbeg, const DenseMatrix &A)
 
void AddSubMatrix (const Array< int > &idx, const DenseMatrix &A)
 (*this)(idx[i],idx[j]) += A(i,j) More...
 
void AddSubMatrix (const Array< int > &idx_i, const Array< int > &idx_j, const DenseMatrix &A)
 (*this)(idx_i[i],idx_j[j]) += A(i,j) More...
 
void AddSubMatrix (int ibeg, const DenseMatrix &A)
 
void AddSubMatrix (int ibeg, int jbeg, const DenseMatrix &A)
 
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
 Prints operator in Matlab format. More...
 
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...
 
std::size_t 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...
 
void Swap (DenseMatrix &other)
 
virtual ~DenseMatrix ()
 Destroys dense matrix. More...
 
virtual void Mult (const Vector &x, Vector &y) const=0
 Operator application: y=A(x). More...
 
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 error. 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 void ArrayMult (const Array< const Vector *> &X, Array< Vector *> &Y) const
 Operator application on a matrix: Y=A(X). More...
 
virtual void ArrayMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y) const
 Action of the transpose operator on a matrix: Y=A^t(X). More...
 
virtual void ArrayAddMult (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const
 Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X). More...
 
virtual void ArrayAddMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const
 Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X). More...
 
virtual OperatorGetGradient (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 void AssembleDiagonal (Vector &diag) const
 Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed. More...
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
virtual const OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. More...
 
virtual const OperatorGetOutputRestriction () 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, 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::Operator
enum  DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP }
 Defines operator diagonal policy upon elimination of rows and/or columns. More...
 
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,
  Complex_DenseMat, MFEM_Block_Matrix, MFEM_Block_Operator
}
 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...
 
OperatorSetupRAP (const Operator *Pi, const Operator *Po)
 Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". 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...
 

Detailed Description

Data type dense matrix using column-major storage.

Definition at line 23 of file densemat.hpp.

Constructor & Destructor Documentation

◆ DenseMatrix() [1/7]

mfem::DenseMatrix::DenseMatrix ( )

Default constructor for DenseMatrix. Sets data = NULL and height = width = 0.

Definition at line 89 of file densemat.cpp.

◆ DenseMatrix() [2/7]

mfem::DenseMatrix::DenseMatrix ( const DenseMatrix m)

Copy constructor.

Definition at line 91 of file densemat.cpp.

◆ DenseMatrix() [3/7]

mfem::DenseMatrix::DenseMatrix ( int  s)
explicit

Creates square matrix of size s.

Definition at line 102 of file densemat.cpp.

◆ DenseMatrix() [4/7]

mfem::DenseMatrix::DenseMatrix ( int  m,
int  n 
)

Creates rectangular matrix of size m x n.

Definition at line 112 of file densemat.cpp.

◆ DenseMatrix() [5/7]

mfem::DenseMatrix::DenseMatrix ( const DenseMatrix mat,
char  ch 
)

Creates rectangular matrix equal to the transpose of mat.

Definition at line 124 of file densemat.cpp.

◆ DenseMatrix() [6/7]

mfem::DenseMatrix::DenseMatrix ( double *  d,
int  h,
int  w 
)
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.

◆ DenseMatrix() [7/7]

template<int M, int N>
mfem::DenseMatrix::DenseMatrix ( const double(&)  values[M][N])
inlineexplicit

Create a dense matrix using a braced initializer list The inner lists correspond to rows of the matrix

Definition at line 64 of file densemat.hpp.

◆ ~DenseMatrix()

mfem::DenseMatrix::~DenseMatrix ( )
virtual

Destroys dense matrix.

Definition at line 2313 of file densemat.cpp.

Member Function Documentation

◆ Add() [1/2]

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 580 of file densemat.cpp.

◆ Add() [2/2]

void mfem::DenseMatrix::Add ( const double  c,
const double *  A 
)

Adds the matrix A multiplied by the number c to the matrix, assuming A has the same dimensions and uses column-major layout.

Definition at line 591 of file densemat.cpp.

◆ AddMatrix() [1/2]

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 1722 of file densemat.cpp.

◆ AddMatrix() [2/2]

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 1752 of file densemat.cpp.

◆ AddMult()

void mfem::DenseMatrix::AddMult ( const Vector x,
Vector y,
const double  a = 1.0 
) const
virtual

y += a * A.x

Reimplemented from mfem::Operator.

Definition at line 251 of file densemat.cpp.

◆ AddMult_a()

void mfem::DenseMatrix::AddMult_a ( double  a,
const Vector x,
Vector y 
) const

y += a * A.x

Definition at line 298 of file densemat.cpp.

◆ AddMultTranspose()

void mfem::DenseMatrix::AddMultTranspose ( const Vector x,
Vector y,
const double  a = 1.0 
) const
virtual

y += a * A^t x

Reimplemented from mfem::Operator.

Definition at line 274 of file densemat.cpp.

◆ AddMultTranspose_a()

void mfem::DenseMatrix::AddMultTranspose_a ( double  a,
const Vector x,
Vector y 
) const

y += a * A^t x

Definition at line 316 of file densemat.cpp.

◆ AddSubMatrix() [1/4]

void mfem::DenseMatrix::AddSubMatrix ( const Array< int > &  idx,
const DenseMatrix A 
)

(*this)(idx[i],idx[j]) += A(i,j)

Definition at line 1999 of file densemat.cpp.

◆ AddSubMatrix() [2/4]

void mfem::DenseMatrix::AddSubMatrix ( const Array< int > &  idx_i,
const Array< int > &  idx_j,
const DenseMatrix A 
)

(*this)(idx_i[i],idx_j[j]) += A(i,j)

Definition at line 2027 of file densemat.cpp.

◆ AddSubMatrix() [3/4]

void mfem::DenseMatrix::AddSubMatrix ( int  ibeg,
const DenseMatrix A 
)

Add the submatrix A to this with row and column offset ibeg

Definition at line 2058 of file densemat.cpp.

◆ AddSubMatrix() [4/4]

void mfem::DenseMatrix::AddSubMatrix ( int  ibeg,
int  jbeg,
const DenseMatrix A 
)

Add the submatrix A to this with row and column offsets ibeg and jbeg respectively

Definition at line 2084 of file densemat.cpp.

◆ AddToVector()

void mfem::DenseMatrix::AddToVector ( int  offset,
Vector v 
) const

Add the matrix 'data' to the Vector 'v' at the given 'offset'.

Definition at line 2112 of file densemat.cpp.

◆ AdjustDofDirection()

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 2134 of file densemat.cpp.

◆ CalcEigenvalues()

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 1292 of file densemat.cpp.

◆ CalcSingularvalue()

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 1267 of file densemat.cpp.

◆ CheckFinite()

int mfem::DenseMatrix::CheckFinite ( ) const
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 452 of file densemat.hpp.

◆ Clear()

void mfem::DenseMatrix::Clear ( )
inline

Delete the matrix data array (if owned) and reset the matrix state.

Definition at line 98 of file densemat.hpp.

◆ ClearExternalData()

void mfem::DenseMatrix::ClearExternalData ( )
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 95 of file densemat.hpp.

◆ CopyCols()

void mfem::DenseMatrix::CopyCols ( const DenseMatrix A,
int  col1,
int  col2 
)

Copy columns col1 through col2 from A to *this.

Definition at line 1578 of file densemat.cpp.

◆ CopyExceptMN()

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 1696 of file densemat.cpp.

◆ CopyMN() [1/3]

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 1591 of file densemat.cpp.

◆ CopyMN() [2/3]

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 1604 of file densemat.cpp.

◆ CopyMN() [3/3]

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 1630 of file densemat.cpp.

◆ CopyMNDiag() [1/2]

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 1661 of file densemat.cpp.

◆ CopyMNDiag() [2/2]

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 1678 of file densemat.cpp.

◆ CopyMNt()

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 1617 of file densemat.cpp.

◆ CopyRows()

void mfem::DenseMatrix::CopyRows ( const DenseMatrix A,
int  row1,
int  row2 
)

Copy rows row1 through row2 from A to *this.

Definition at line 1565 of file densemat.cpp.

◆ Data()

double* mfem::DenseMatrix::Data ( ) const
inline

Returns the matrix data array.

Definition at line 111 of file densemat.hpp.

◆ Det()

double mfem::DenseMatrix::Det ( ) const

Calculates the determinant of the matrix (optimized for 2x2, 3x3, and 4x4 matrices)

Definition at line 487 of file densemat.cpp.

◆ Diag() [1/2]

void mfem::DenseMatrix::Diag ( double  c,
int  n 
)

Creates n x n diagonal matrix with diagonal elements c.

Definition at line 1389 of file densemat.cpp.

◆ Diag() [2/2]

void mfem::DenseMatrix::Diag ( double *  diag,
int  n 
)

Creates n x n diagonal matrix with diagonal given by diag.

Definition at line 1404 of file densemat.cpp.

◆ Eigensystem() [1/2]

void mfem::DenseMatrix::Eigensystem ( Vector ev,
DenseMatrix evect 
)
inline

Compute eigenvalues and eigenvectors of A x = ev x where A = *this.

Definition at line 280 of file densemat.hpp.

◆ Eigensystem() [2/2]

void mfem::DenseMatrix::Eigensystem ( DenseMatrix b,
Vector ev,
DenseMatrix evect 
)
inline

Compute generalized eigenvalues and eigenvectors of A x = ev B x, where A = *this

Definition at line 294 of file densemat.hpp.

◆ Eigenvalues() [1/4]

void mfem::DenseMatrix::Eigenvalues ( Vector ev)
inline

Compute eigenvalues of A x = ev x where A = *this.

Definition at line 272 of file densemat.hpp.

◆ Eigenvalues() [2/4]

void mfem::DenseMatrix::Eigenvalues ( Vector ev,
DenseMatrix evect 
)
inline

Compute eigenvalues and eigenvectors of A x = ev x where A = *this.

Definition at line 276 of file densemat.hpp.

◆ Eigenvalues() [3/4]

void mfem::DenseMatrix::Eigenvalues ( DenseMatrix b,
Vector ev 
)
inline

Compute generalized eigenvalues and eigenvectors of A x = ev B x, where A = *this

Definition at line 285 of file densemat.hpp.

◆ Eigenvalues() [4/4]

void mfem::DenseMatrix::Eigenvalues ( DenseMatrix b,
Vector ev,
DenseMatrix evect 
)
inline

Compute generalized eigenvalues of A x = ev B x, where A = *this.

Definition at line 289 of file densemat.hpp.

◆ Elem() [1/2]

double & mfem::DenseMatrix::Elem ( int  i,
int  j 
)
virtual

Returns reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 162 of file densemat.cpp.

◆ Elem() [2/2]

const double & mfem::DenseMatrix::Elem ( int  i,
int  j 
) const
virtual

Returns constant reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 167 of file densemat.cpp.

◆ FNorm()

double mfem::DenseMatrix::FNorm ( ) const
inline

Compute the Frobenius norm of the matrix.

Definition at line 266 of file densemat.hpp.

◆ FNorm2()

double mfem::DenseMatrix::FNorm2 ( ) const
inline

Compute the square of the Frobenius norm of the matrix.

Definition at line 269 of file densemat.hpp.

◆ GetColumn() [1/3]

void mfem::DenseMatrix::GetColumn ( int  c,
Vector col 
) const

Definition at line 1330 of file densemat.cpp.

◆ GetColumn() [2/3]

double* mfem::DenseMatrix::GetColumn ( int  col)
inline

Definition at line 309 of file densemat.hpp.

◆ GetColumn() [3/3]

const double* mfem::DenseMatrix::GetColumn ( int  col) const
inline

Definition at line 310 of file densemat.hpp.

◆ GetColumnReference()

void mfem::DenseMatrix::GetColumnReference ( int  c,
Vector col 
)
inline

Definition at line 312 of file densemat.hpp.

◆ GetData()

double* mfem::DenseMatrix::GetData ( ) const
inline

Returns the matrix data array.

Definition at line 115 of file densemat.hpp.

◆ GetDiag()

void mfem::DenseMatrix::GetDiag ( Vector d) const

Returns the diagonal of the matrix.

Definition at line 1344 of file densemat.cpp.

◆ GetFromVector()

void mfem::DenseMatrix::GetFromVector ( int  offset,
const Vector v 
)

Get the matrix 'data' from the Vector 'v' at the given 'offset'.

Definition at line 2123 of file densemat.cpp.

◆ Getl1Diag()

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 1358 of file densemat.cpp.

◆ GetMemory() [1/2]

Memory<double>& mfem::DenseMatrix::GetMemory ( )
inline

Definition at line 117 of file densemat.hpp.

◆ GetMemory() [2/2]

const Memory<double>& mfem::DenseMatrix::GetMemory ( ) const
inline

Definition at line 118 of file densemat.hpp.

◆ GetRow()

void mfem::DenseMatrix::GetRow ( int  r,
Vector row 
) const

Definition at line 1314 of file densemat.cpp.

◆ GetRowSums()

void mfem::DenseMatrix::GetRowSums ( Vector l) const

Compute the row sums of the DenseMatrix.

Definition at line 1375 of file densemat.cpp.

◆ GetSubMatrix() [1/4]

void mfem::DenseMatrix::GetSubMatrix ( const Array< int > &  idx,
DenseMatrix A 
) const

Get the square submatrix which corresponds to the given indices idx. Note: the A matrix will be resized to accommodate the data

Definition at line 1782 of file densemat.cpp.

◆ GetSubMatrix() [2/4]

void mfem::DenseMatrix::GetSubMatrix ( const Array< int > &  idx_i,
const Array< int > &  idx_j,
DenseMatrix A 
) const

Get the rectangular submatrix which corresponds to the given indices idx_i and idx_j. Note: the A matrix will be resized to accommodate the data

Definition at line 1803 of file densemat.cpp.

◆ GetSubMatrix() [3/4]

void mfem::DenseMatrix::GetSubMatrix ( int  ibeg,
int  iend,
DenseMatrix A 
)

Get the square submatrix which corresponds to the range [ ibeg, iend ). Note: the A matrix will be resized to accommodate the data

Definition at line 1829 of file densemat.cpp.

◆ GetSubMatrix() [4/4]

void mfem::DenseMatrix::GetSubMatrix ( int  ibeg,
int  iend,
int  jbeg,
int  jend,
DenseMatrix A 
)

Get the square submatrix which corresponds to the range i ∈ [ ibeg, iend ) and j ∈ [ jbeg, jend ) Note: the A matrix will be resized to accommodate the data

Definition at line 1853 of file densemat.cpp.

◆ GradToCurl()

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. In 2D this computes the scalar-valued curl of a 2D vector field

Definition at line 1477 of file densemat.cpp.

◆ GradToDiv()

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 1550 of file densemat.cpp.

◆ GradToVectorCurl2D()

void mfem::DenseMatrix::GradToVectorCurl2D ( DenseMatrix curl)

Given a DShape matrix (from a scalar FE), stored in *this, returns the CurlShape matrix. This computes the vector-valued curl of a scalar field. this is N by 2 matrix and curl is N by 2 matrix as well.

Definition at line 1536 of file densemat.cpp.

◆ HostRead()

const double* mfem::DenseMatrix::HostRead ( ) const
inline

Shortcut for mfem::Read(GetMemory(), TotalSize(), false).

Definition at line 470 of file densemat.hpp.

◆ HostReadWrite()

double* mfem::DenseMatrix::HostReadWrite ( )
inline

Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false).

Definition at line 486 of file densemat.hpp.

◆ HostWrite()

double* mfem::DenseMatrix::HostWrite ( )
inline

Shortcut for mfem::Write(GetMemory(), TotalSize(), false).

Definition at line 478 of file densemat.hpp.

◆ InnerProduct() [1/2]

double mfem::DenseMatrix::InnerProduct ( const double *  x,
const double *  y 
) const

Compute y^t A x.

Definition at line 335 of file densemat.cpp.

◆ InnerProduct() [2/2]

double mfem::DenseMatrix::InnerProduct ( const Vector x,
const Vector y 
) const
inline

Compute y^t A x.

Definition at line 198 of file densemat.hpp.

◆ Inverse()

MatrixInverse * mfem::DenseMatrix::Inverse ( ) const
virtual

Returns a pointer to the inverse matrix.

Implements mfem::Matrix.

Definition at line 482 of file densemat.cpp.

◆ Invert()

void mfem::DenseMatrix::Invert ( )

Replaces the current matrix with its inverse.

Definition at line 678 of file densemat.cpp.

◆ InvLeftScaling()

void mfem::DenseMatrix::InvLeftScaling ( const Vector s)

InvLeftScaling this = diag(1./s) * this.

Definition at line 366 of file densemat.cpp.

◆ InvRightScaling()

void mfem::DenseMatrix::InvRightScaling ( const Vector s)

InvRightScaling: this = this * diag(1./s);.

Definition at line 394 of file densemat.cpp.

◆ InvSymmetricScaling()

void mfem::DenseMatrix::InvSymmetricScaling ( const Vector s)

InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))

Definition at line 436 of file densemat.cpp.

◆ LeftScaling()

void mfem::DenseMatrix::LeftScaling ( const Vector s)

LeftScaling this = diag(s) * this.

Definition at line 353 of file densemat.cpp.

◆ Lump()

void mfem::DenseMatrix::Lump ( )

Definition at line 1463 of file densemat.cpp.

◆ MaxMaxNorm()

double mfem::DenseMatrix::MaxMaxNorm ( ) const

Compute the norm ||A|| = max_{ij} |A_{ij}|.

Definition at line 846 of file densemat.cpp.

◆ MemoryUsage()

std::size_t mfem::DenseMatrix::MemoryUsage ( ) const
inline

Definition at line 463 of file densemat.hpp.

◆ Mult() [1/5]

void mfem::DenseMatrix::Mult ( const double *  x,
double *  y 
) const

Matrix vector multiplication.

Definition at line 172 of file densemat.cpp.

◆ Mult() [2/5]

void mfem::DenseMatrix::Mult ( const double *  x,
Vector y 
) const

Matrix vector multiplication.

Definition at line 177 of file densemat.cpp.

◆ Mult() [3/5]

void mfem::DenseMatrix::Mult ( const Vector x,
double *  y 
) const

Matrix vector multiplication.

Definition at line 184 of file densemat.cpp.

◆ Mult() [4/5]

void mfem::DenseMatrix::Mult ( const Vector x,
Vector y 
) const
virtual

Matrix vector multiplication.

Implements mfem::Operator.

Definition at line 191 of file densemat.cpp.

◆ Mult() [5/5]

virtual void mfem::Operator::Mult

Operator application: y=A(x).

◆ MultTranspose() [1/5]

void mfem::DenseMatrix::MultTranspose ( const double *  x,
double *  y 
) const

Multiply a vector with the transpose matrix.

Definition at line 214 of file densemat.cpp.

◆ MultTranspose() [2/5]

void mfem::DenseMatrix::MultTranspose ( const double *  x,
Vector y 
) const

Multiply a vector with the transpose matrix.

Definition at line 229 of file densemat.cpp.

◆ MultTranspose() [3/5]

void mfem::DenseMatrix::MultTranspose ( const Vector x,
double *  y 
) const

Multiply a vector with the transpose matrix.

Definition at line 236 of file densemat.cpp.

◆ MultTranspose() [4/5]

void mfem::DenseMatrix::MultTranspose ( const Vector x,
Vector y 
) const
virtual

Multiply a vector with the transpose matrix.

Reimplemented from mfem::Operator.

Definition at line 243 of file densemat.cpp.

◆ MultTranspose() [5/5]

virtual void mfem::Operator::MultTranspose
inline

Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an error.

Definition at line 93 of file operator.hpp.

◆ Neg()

void mfem::DenseMatrix::Neg ( )

(*this) = -(*this)

Definition at line 669 of file densemat.cpp.

◆ Norm2() [1/2]

void mfem::DenseMatrix::Norm2 ( double *  v) const

Take the 2-norm of the columns of A and store in v.

Definition at line 833 of file densemat.cpp.

◆ Norm2() [2/2]

void mfem::DenseMatrix::Norm2 ( Vector v) const
inline

Take the 2-norm of the columns of A and store in v.

Definition at line 256 of file densemat.hpp.

◆ operator()() [1/2]

double & mfem::DenseMatrix::operator() ( int  i,
int  j 
)
inline

Returns reference to a_{ij}.

Definition at line 1288 of file densemat.hpp.

◆ operator()() [2/2]

const double & mfem::DenseMatrix::operator() ( int  i,
int  j 
) const
inline

Returns constant reference to a_{ij}.

Definition at line 1294 of file densemat.hpp.

◆ operator*()

double mfem::DenseMatrix::operator* ( const DenseMatrix m) const

Matrix inner product: tr(A^t B)

Definition at line 199 of file densemat.cpp.

◆ operator*=()

DenseMatrix & mfem::DenseMatrix::operator*= ( double  c)

Definition at line 659 of file densemat.cpp.

◆ operator+=() [1/2]

DenseMatrix & mfem::DenseMatrix::operator+= ( const double *  m)

Definition at line 633 of file densemat.cpp.

◆ operator+=() [2/2]

DenseMatrix & mfem::DenseMatrix::operator+= ( const DenseMatrix m)

Definition at line 639 of file densemat.cpp.

◆ operator-=()

DenseMatrix & mfem::DenseMatrix::operator-= ( const DenseMatrix m)

Definition at line 646 of file densemat.cpp.

◆ operator=() [1/3]

DenseMatrix & mfem::DenseMatrix::operator= ( double  c)

Sets the matrix elements equal to constant c.

Definition at line 600 of file densemat.cpp.

◆ operator=() [2/3]

DenseMatrix & mfem::DenseMatrix::operator= ( const double *  d)

Copy the matrix entries from the given array.

Definition at line 610 of file densemat.cpp.

◆ operator=() [3/3]

DenseMatrix & mfem::DenseMatrix::operator= ( const DenseMatrix m)

Sets the matrix size and elements equal to those of m.

Definition at line 620 of file densemat.cpp.

◆ OwnsData()

bool mfem::DenseMatrix::OwnsData ( ) const
inline

Return the DenseMatrix data (host pointer) ownership flag.

Definition at line 121 of file densemat.hpp.

◆ Print()

void mfem::DenseMatrix::Print ( std::ostream &  out = mfem::out,
int  width_ = 4 
) const
virtual

Prints matrix to stream out.

Reimplemented from mfem::Matrix.

Definition at line 2221 of file densemat.cpp.

◆ PrintMatlab()

void mfem::DenseMatrix::PrintMatlab ( std::ostream &  out = mfem::out) const
virtual

Prints operator in Matlab format.

Reimplemented from mfem::Operator.

Definition at line 2247 of file densemat.cpp.

◆ PrintT()

void mfem::DenseMatrix::PrintT ( std::ostream &  out = mfem::out,
int  width_ = 4 
) const
virtual

Prints the transpose matrix to stream out.

Definition at line 2266 of file densemat.cpp.

◆ Rank()

int mfem::DenseMatrix::Rank ( double  tol) const

Definition at line 1252 of file densemat.cpp.

◆ Read()

const double* mfem::DenseMatrix::Read ( bool  on_dev = true) const
inline

Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).

Definition at line 466 of file densemat.hpp.

◆ ReadWrite()

double* mfem::DenseMatrix::ReadWrite ( bool  on_dev = true)
inline

Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev).

Definition at line 482 of file densemat.hpp.

◆ Reset()

void mfem::DenseMatrix::Reset ( double *  d,
int  h,
int  w 
)
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 90 of file densemat.hpp.

◆ RightScaling()

void mfem::DenseMatrix::RightScaling ( const Vector s)

RightScaling: this = this * diag(s);.

Definition at line 379 of file densemat.cpp.

◆ Set() [1/2]

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 571 of file densemat.cpp.

◆ Set() [2/2]

void mfem::DenseMatrix::Set ( double  alpha,
const DenseMatrix A 
)
inline

Set the matrix to alpha * A.

Definition at line 220 of file densemat.hpp.

◆ SetCol() [1/3]

void mfem::DenseMatrix::SetCol ( int  c,
const double *  col 
)

Definition at line 2192 of file densemat.cpp.

◆ SetCol() [2/3]

void mfem::DenseMatrix::SetCol ( int  c,
const Vector col 
)

Definition at line 2201 of file densemat.cpp.

◆ SetCol() [3/3]

void mfem::DenseMatrix::SetCol ( int  col,
double  value 
)

Set all entries of a column to the specified value.

Definition at line 2169 of file densemat.cpp.

◆ SetRow() [1/3]

void mfem::DenseMatrix::SetRow ( int  r,
const double *  row 
)

Definition at line 2177 of file densemat.cpp.

◆ SetRow() [2/3]

void mfem::DenseMatrix::SetRow ( int  r,
const Vector row 
)

Definition at line 2186 of file densemat.cpp.

◆ SetRow() [3/3]

void mfem::DenseMatrix::SetRow ( int  row,
double  value 
)

Set all entries of a row to the specified value.

Definition at line 2161 of file densemat.cpp.

◆ SetSize() [1/2]

void mfem::DenseMatrix::SetSize ( int  s)
inline

Change the size of the DenseMatrix to s x s.

Definition at line 105 of file densemat.hpp.

◆ SetSize() [2/2]

void mfem::DenseMatrix::SetSize ( int  h,
int  w 
)

Change the size of the DenseMatrix to h x w.

Definition at line 143 of file densemat.cpp.

◆ SetSubMatrix() [1/4]

void mfem::DenseMatrix::SetSubMatrix ( const Array< int > &  idx,
const DenseMatrix A 
)

Set (*this)(idx[i],idx[j]) = A(i,j)

Definition at line 1886 of file densemat.cpp.

◆ SetSubMatrix() [2/4]

void mfem::DenseMatrix::SetSubMatrix ( const Array< int > &  idx_i,
const Array< int > &  idx_j,
const DenseMatrix A 
)

Set (*this)(idx_i[i],idx_j[j]) = A(i,j)

Definition at line 1915 of file densemat.cpp.

◆ SetSubMatrix() [3/4]

void mfem::DenseMatrix::SetSubMatrix ( int  ibeg,
const DenseMatrix A 
)

Set a submatrix of (this) to the given matrix A with row and column offset ibeg

Definition at line 1945 of file densemat.cpp.

◆ SetSubMatrix() [4/4]

void mfem::DenseMatrix::SetSubMatrix ( int  ibeg,
int  jbeg,
const DenseMatrix A 
)

Set a submatrix of (this) to the given matrix A with row and column offset ibeg and jbeg respectively

Definition at line 1971 of file densemat.cpp.

◆ SingularValues()

void mfem::DenseMatrix::SingularValues ( Vector sv) const

Definition at line 1212 of file densemat.cpp.

◆ Size()

int mfem::DenseMatrix::Size ( ) const
inline

For backward compatibility define Size to be synonym of Width()

Definition at line 102 of file densemat.hpp.

◆ SquareRootInverse()

void mfem::DenseMatrix::SquareRootInverse ( )

Replaces the current matrix with its square root inverse.

Definition at line 788 of file densemat.cpp.

◆ Swap()

void mfem::DenseMatrix::Swap ( DenseMatrix other)

Definition at line 2306 of file densemat.cpp.

◆ SymmetricScaling()

void mfem::DenseMatrix::SymmetricScaling ( const Vector s)

SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))

Definition at line 408 of file densemat.cpp.

◆ Symmetrize()

void mfem::DenseMatrix::Symmetrize ( )

(*this) = 1/2 ((*this) + (*this)^t)

Definition at line 1452 of file densemat.cpp.

◆ TestInversion()

void mfem::DenseMatrix::TestInversion ( )

Invert and print the numerical conditioning of the inversion.

Definition at line 2292 of file densemat.cpp.

◆ Threshold()

void mfem::DenseMatrix::Threshold ( double  eps)

Replace small entries, abs(a_ij) <= eps, with zero.

Definition at line 2207 of file densemat.cpp.

◆ Trace()

double mfem::DenseMatrix::Trace ( ) const

Trace of a square matrix.

Definition at line 463 of file densemat.cpp.

◆ Transpose() [1/2]

void mfem::DenseMatrix::Transpose ( )

(*this) = (*this)^t

Definition at line 1419 of file densemat.cpp.

◆ Transpose() [2/2]

void mfem::DenseMatrix::Transpose ( const DenseMatrix A)

(*this) = A^t

Definition at line 1441 of file densemat.cpp.

◆ UseExternalData()

void mfem::DenseMatrix::UseExternalData ( double *  d,
int  h,
int  w 
)
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 80 of file densemat.hpp.

◆ Weight()

double mfem::DenseMatrix::Weight ( ) const

Definition at line 544 of file densemat.cpp.

◆ Write()

double* mfem::DenseMatrix::Write ( bool  on_dev = true)
inline

Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).

Definition at line 474 of file densemat.hpp.

Friends And Related Function Documentation

◆ DenseMatrixInverse

friend class DenseMatrixInverse
friend

Definition at line 26 of file densemat.hpp.

◆ DenseTensor

friend class DenseTensor
friend

Definition at line 25 of file densemat.hpp.


The documentation for this class was generated from the following files: