|
MFEM v2.0
|
Data type dense matrix. More...
#include <densemat.hpp>


Public Member Functions | |
| DenseMatrix () | |
| DenseMatrix (const DenseMatrix &) | |
| Copy constructor. | |
| DenseMatrix (int s) | |
| Creates square matrix of size s. | |
| DenseMatrix (int m, int n) | |
| Creates rectangular matrix of size n and height m. | |
| DenseMatrix (const DenseMatrix &mat, char ch) | |
| Creates rectangular matrix equal to the transpose of mat. | |
| DenseMatrix (double *d, int h, int w) | |
| void | UseExternalData (double *d, int h, int w) |
| void | ClearExternalData () |
| void | SetSize (int s) |
| If the matrix is not a square matrix of size s then recreate it. | |
| void | SetSize (int h, int w) |
| If the matrix is not a matrix of size (h x w) then recreate it. | |
| int | Height () const |
| Returns the height of the matrix. | |
| int | Width () const |
| Returns the width of the matrix. | |
| double * | Data () const |
| Returns vector of the elements. | |
| double & | operator() (int i, int j) |
| Returns reference to a_{ij}. Index i, j = 0 .. size-1. | |
| const double & | operator() (int i, int j) const |
| Returns constant reference to a_{ij}. Index i, j = 0 .. size-1. | |
| double | operator* (const DenseMatrix &m) const |
| Matrix inner product: tr(A^t B) | |
| double | Trace () const |
| Trace of a square matrix. | |
| virtual double & | Elem (int i, int j) |
| Returns reference to a_{ij}. Index i, j = 0 .. size-1. | |
| virtual const double & | Elem (int i, int j) const |
| Returns constant reference to a_{ij}. Index i, j = 0 .. size-1. | |
| void | Mult (const double *x, double *y) const |
| Matrix vector multiplication. | |
| virtual void | Mult (const Vector &x, Vector &y) const |
| Matrix vector multiplication. | |
| virtual void | MultTranspose (const Vector &x, Vector &y) const |
| Multiply a vector with the transpose matrix. | |
| void | AddMult (const Vector &x, Vector &y) const |
| y += A.x | |
| double | InnerProduct (const double *x, const double *y) const |
| Compute y^t A x. | |
| double | InnerProduct (const Vector &x, const Vector &y) const |
| Compute y^t A x. | |
| virtual MatrixInverse * | Inverse () const |
| Returns a pointer to the inverse matrix. | |
| void | Invert () |
| Replaces the current matrix with its inverse. | |
| double | Det () const |
| Calculates the determinant of the matrix (for 2x2 or 3x3 matrices) | |
| double | Weight () const |
| void | Add (const double c, DenseMatrix &A) |
| Adds the matrix A multiplied by the number c to the matrix. | |
| DenseMatrix & | operator= (double c) |
| Sets the matrix elements equal to constant c. | |
| DenseMatrix & | operator= (const DenseMatrix &m) |
| Sets the matrix size and elements equal to those of m. | |
| DenseMatrix & | operator+= (DenseMatrix &m) |
| DenseMatrix & | operator-= (DenseMatrix &m) |
| DenseMatrix & | operator*= (double c) |
| void | Neg () |
| (*this) = -(*this) | |
| void | Norm2 (double *v) const |
| Take the 2-norm of the columns of A and store in v. | |
| double | FNorm () const |
| void | Eigenvalues (Vector &ev) |
| void | Eigenvalues (Vector &ev, DenseMatrix &evect) |
| void | Eigensystem (Vector &ev, DenseMatrix &evect) |
| void | SingularValues (Vector &sv) const |
| double | CalcSingularvalue (const int i) const |
| Return the i-th singular value (decreasing order) of a 2x2 or 3x3 matrix. | |
| void | CalcEigenvalues (double *lambda, double *vec) const |
| void | GetColumn (int c, Vector &col) |
| void | GetColumnReference (int c, Vector &col) |
| void | Diag (double c, int n) |
| Creates n x n diagonal matrix with diagonal elements c. | |
| void | Diag (double *diag, int n) |
| Creates n x n diagonal matrix with diagonal given by diag. | |
| void | Transpose () |
| (*this) = (*this)^t | |
| void | Transpose (DenseMatrix &A) |
| (*this) = A^t | |
| void | Symmetrize () |
| (*this) = 1/2 ((*this) + (*this)^t) | |
| void | Lump () |
| void | GradToCurl (DenseMatrix &curl) |
| void | GradToDiv (Vector &div) |
| void | CopyRows (DenseMatrix &A, int row1, int row2) |
| Copy rows row1 through row2 from A to *this. | |
| void | CopyCols (DenseMatrix &A, int col1, int col2) |
| Copy columns col1 through col2 from A to *this. | |
| void | CopyMN (DenseMatrix &A, int m, int n, int Aro, int Aco) |
| Copy the m x n submatrix of A at (Aro)ffset, (Aco)loffset to *this. | |
| void | CopyMN (DenseMatrix &A, int row_offset, int col_offset) |
| Copy matrix A to the location in *this at row_offset, col_offset. | |
| void | CopyMNt (DenseMatrix &A, int row_offset, int col_offset) |
| Copy matrix A^t to the location in *this at row_offset, 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. | |
| 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. | |
| 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. | |
| void | AddMatrix (double a, DenseMatrix &A, int ro, int co) |
| Perform (ro+i,co+j)+=a*A(i,j) for 0<=i<A.Height, 0<=j<A.Width. | |
| void | AddToVector (int offset, Vector &v) const |
| Add the matrix 'data' to the Vector 'v' at the given 'offset'. | |
| void | GetFromVector (int offset, const Vector &v) |
| Get the matrix 'data' from the Vector 'v' at the given 'offset'. | |
| void | AdjustDofDirection (Array< int > &dofs) |
| virtual void | Print (ostream &out=cout, int width=4) const |
| Prints matrix to stream out. | |
| virtual void | PrintT (ostream &out=cout, int width=4) const |
| Prints the transpose matrix to stream out. | |
| void | TestInversion () |
| Invert and print the numerical conditioning of the inversion. | |
| virtual | ~DenseMatrix () |
| Destroys dense matrix. | |
Private Member Functions | |
| void | Eigensystem (Vector &ev, DenseMatrix *evect=NULL) |
Private Attributes | |
| int | height |
| double * | data |
Friends | |
| class | DenseTensor |
| class | DenseMatrixInverse |
| void | Mult (const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a) |
| Matrix matrix multiplication. A = B * C. | |
Data type dense matrix.
Definition at line 16 of file densemat.hpp.
| DenseMatrix::DenseMatrix | ( | ) |
Default constructor for DenseMatrix. Sets data = NULL size = height = 0
Definition at line 26 of file densemat.cpp.
| DenseMatrix::DenseMatrix | ( | const DenseMatrix & | m | ) |
Copy constructor.
Definition at line 32 of file densemat.cpp.
References data, height, and Operator::size.
| DenseMatrix::DenseMatrix | ( | int | s | ) | [explicit] |
Creates square matrix of size s.
Definition at line 41 of file densemat.cpp.
| DenseMatrix::DenseMatrix | ( | int | m, |
| int | n | ||
| ) |
Creates rectangular matrix of size n and height m.
Definition at line 52 of file densemat.cpp.
| DenseMatrix::DenseMatrix | ( | const DenseMatrix & | mat, |
| char | ch | ||
| ) |
Creates rectangular matrix equal to the transpose of mat.
Definition at line 63 of file densemat.cpp.
References data, height, and Operator::size.
| DenseMatrix::DenseMatrix | ( | double * | d, |
| int | h, | ||
| int | w | ||
| ) | [inline] |
Definition at line 48 of file densemat.hpp.
| DenseMatrix::~DenseMatrix | ( | ) | [virtual] |
| void DenseMatrix::Add | ( | const double | c, |
| DenseMatrix & | A | ||
| ) |
Adds the matrix A multiplied by the number c to the matrix.
Definition at line 273 of file densemat.cpp.
References Height(), and Operator::Size().
| void 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 1637 of file densemat.cpp.
References data, Height(), mfem_error(), and Width().
Referenced by VectorMassIntegrator::AssembleElementMatrix().
| void DenseMatrix::AddMatrix | ( | double | a, |
| 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 1663 of file densemat.cpp.
References data, Height(), mfem_error(), and Width().
y += A.x
Definition at line 176 of file densemat.cpp.
References height, mfem_error(), Vector::Size(), and Operator::size.
Referenced by VectorFEDomainLFIntegrator::AssembleRHSElementVect().
| void DenseMatrix::AddToVector | ( | int | offset, |
| Vector & | v | ||
| ) | const |
Add the matrix 'data' to the Vector 'v' at the given 'offset'.
Definition at line 1689 of file densemat.cpp.
References data, Vector::GetData(), height, and Operator::size.
| void 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 1707 of file densemat.cpp.
References Height(), mfem_error(), Array< T >::Size(), and Width().
| void 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 1120 of file densemat.cpp.
References data, Height(), mfem_error(), and Width().
| double DenseMatrix::CalcSingularvalue | ( | const int | i | ) | const |
Return the i-th singular value (decreasing order) of a 2x2 or 3x3 matrix.
Definition at line 738 of file densemat.cpp.
References data, Height(), mfem_error(), and Width().
Referenced by Mesh::GetElementSize(), and Mesh::PrintCharacteristics().
| void DenseMatrix::ClearExternalData | ( | ) | [inline] |
Definition at line 52 of file densemat.hpp.
References data, height, and Operator::size.
Referenced by BilinearForm::ComputeElementMatrices(), and DenseTensor::~DenseTensor().
| void DenseMatrix::CopyCols | ( | DenseMatrix & | A, |
| int | col1, | ||
| int | col2 | ||
| ) |
Copy columns col1 through col2 from A to *this.
Definition at line 1570 of file densemat.cpp.
| void DenseMatrix::CopyMN | ( | DenseMatrix & | A, |
| int | m, | ||
| int | n, | ||
| int | Aro, | ||
| int | Aco | ||
| ) |
Copy the m x n submatrix of A at (Aro)ffset, (Aco)loffset to *this.
Definition at line 1579 of file densemat.cpp.
References SetSize().
| void DenseMatrix::CopyMN | ( | DenseMatrix & | A, |
| int | row_offset, | ||
| int | col_offset | ||
| ) |
Copy matrix A to the location in *this at row_offset, col_offset.
Definition at line 1590 of file densemat.cpp.
References data, Height(), and Operator::Size().
| void 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 1610 of file densemat.cpp.
| void 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 1623 of file densemat.cpp.
| void DenseMatrix::CopyMNt | ( | 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 1600 of file densemat.cpp.
References data, Height(), and Operator::Size().
| void DenseMatrix::CopyRows | ( | DenseMatrix & | A, |
| int | row1, | ||
| int | row2 | ||
| ) |
Copy rows row1 through row2 from A to *this.
Definition at line 1561 of file densemat.cpp.
References SetSize(), and Operator::Size().
| double* DenseMatrix::Data | ( | ) | const [inline] |
Returns vector of the elements.
Definition at line 67 of file densemat.hpp.
References data.
Referenced by NURBS1DFiniteElement::CalcDShape(), L2_SegmentElement::CalcDShape(), CalcInverse(), DenseMatrixEigensystem::DenseMatrixEigensystem(), DetOfLinComb(), dsyev_Eigensystem(), dsyevr_Eigensystem(), DenseMatrixEigensystem::Eigenvector(), DenseMatrixSVD::Eval(), DenseMatrixEigensystem::Eval(), MultABt(), and MultAtB().
| double DenseMatrix::Det | ( | ) | const |
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
Definition at line 223 of file densemat.cpp.
References data, Height(), mfem_error(), and Width().
Referenced by CalcInverse(), CalcInverseTranspose(), Mesh::CheckElementOrientation(), Mesh::GetElementSize(), Mesh::PrintCharacteristics(), and Weight().
| void DenseMatrix::Diag | ( | double | c, |
| int | n | ||
| ) |
Creates n x n diagonal matrix with diagonal elements c.
Definition at line 1426 of file densemat.cpp.
| void DenseMatrix::Diag | ( | double * | diag, |
| int | n | ||
| ) |
Creates n x n diagonal matrix with diagonal given by diag.
Definition at line 1437 of file densemat.cpp.
| void DenseMatrix::Eigensystem | ( | Vector & | ev, |
| DenseMatrix & | evect | ||
| ) | [inline] |
| void DenseMatrix::Eigensystem | ( | Vector & | ev, |
| DenseMatrix * | evect = NULL |
||
| ) | [private] |
Definition at line 684 of file densemat.cpp.
References dsyev_Eigensystem(), and mfem_error().
Referenced by Eigenvalues().
| void DenseMatrix::Eigenvalues | ( | Vector & | ev | ) | [inline] |
Definition at line 140 of file densemat.hpp.
References Eigensystem().
| void DenseMatrix::Eigenvalues | ( | Vector & | ev, |
| DenseMatrix & | evect | ||
| ) | [inline] |
Definition at line 143 of file densemat.hpp.
References Eigensystem().
| double & DenseMatrix::Elem | ( | int | i, |
| int | j | ||
| ) | [virtual] |
Returns reference to a_{ij}. Index i, j = 0 .. size-1.
Implements Matrix.
Definition at line 114 of file densemat.cpp.
| const double & DenseMatrix::Elem | ( | int | i, |
| int | j | ||
| ) | const [virtual] |
Returns constant reference to a_{ij}. Index i, j = 0 .. size-1.
Implements Matrix.
Definition at line 119 of file densemat.cpp.
| double DenseMatrix::FNorm | ( | ) | const |
Definition at line 448 of file densemat.cpp.
References data, Height(), and Operator::Size().
Referenced by CalcInverse(), and TestInversion().
| void DenseMatrix::GetColumn | ( | int | c, |
| Vector & | col | ||
| ) |
Definition at line 1412 of file densemat.cpp.
References data, Vector::GetData(), Height(), and Vector::SetSize().
| void DenseMatrix::GetColumnReference | ( | int | c, |
| Vector & | col | ||
| ) | [inline] |
Definition at line 160 of file densemat.hpp.
References data, height, and Vector::SetDataAndSize().
Referenced by VectorCoefficient::Eval(), GridFunction::GetGradients(), and GridFunction::GetVectorValues().
| void DenseMatrix::GetFromVector | ( | int | offset, |
| const Vector & | v | ||
| ) |
Get the matrix 'data' from the Vector 'v' at the given 'offset'.
Definition at line 1698 of file densemat.cpp.
References data, Vector::GetData(), height, and Operator::size.
| void 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 3 matrix, then curl is a 3*N by 3 matrix. The size of curl must be set outside.
Definition at line 1508 of file densemat.cpp.
References Height(), mfem_error(), and Width().
| void 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 1544 of file densemat.cpp.
References data, Vector::GetData(), height, Height(), mfem_error(), Operator::size, Vector::Size(), and Width().
Referenced by ElasticityIntegrator::AssembleElementMatrix(), and VectorDivergenceIntegrator::AssembleElementMatrix2().
| int DenseMatrix::Height | ( | ) | const [inline] |
Returns the height of the matrix.
Definition at line 61 of file densemat.hpp.
References height.
Referenced by Add(), Add(), AddMatrix(), AddMult_a_AAt(), AddMult_a_VVt(), AddMultABt(), AddMultADAt(), AddMultVWt(), AdjustDofDirection(), CalcAdjugate(), CalcAdjugateTranspose(), Linear3DFiniteElement::CalcDShape(), CalcEigenvalues(), CalcInverse(), CalcInverseTranspose(), CalcSingularvalue(), VectorFiniteElement::CalcVShape_ND(), CopyCols(), CopyMN(), CopyMNt(), DenseMatrixSVD::DenseMatrixSVD(), Det(), dsyev_Eigensystem(), dsyevr_Eigensystem(), DenseMatrixEigensystem::Eigenvector(), DenseMatrixSVD::Eval(), FNorm(), GetColumn(), DenseTensor::GetData(), GradToCurl(), GradToDiv(), Invert(), IsoparametricTransformation::Jacobian(), Lump(), Mult_a_AAt(), MultAAt(), MultABt(), MultADAt(), MultAtB(), MultVWt(), Neg(), Norm2(), operator*(), operator*=(), operator=(), Mesh::PrintVTK(), GridFunction::SaveVTK(), SetSize(), SingularValues(), DenseTensor::SizeI(), Symmetrize(), Trace(), IntegrationPointTransformation::Transform(), IsoparametricTransformation::Transform(), Transpose(), and Weight().
| double DenseMatrix::InnerProduct | ( | const double * | x, |
| const double * | y | ||
| ) | const |
Compute y^t A x.
Definition at line 188 of file densemat.cpp.
References height, and Operator::size.
Referenced by VectorFiniteElement::Project_ND(), and VectorFiniteElement::Project_RT().
Compute y^t A x.
Definition at line 103 of file densemat.hpp.
References InnerProduct().
Referenced by InnerProduct().
| MatrixInverse * DenseMatrix::Inverse | ( | ) | const [virtual] |
Returns a pointer to the inverse matrix.
Implements Matrix.
Definition at line 218 of file densemat.cpp.
References DenseMatrixInverse.
| void DenseMatrix::Invert | ( | ) |
Replaces the current matrix with its inverse.
Definition at line 351 of file densemat.cpp.
References data, Height(), mfem_error(), Operator::size, and Operator::Size().
Referenced by Poly_1D::Basis::Basis(), GaussQuad2DFiniteElement::GaussQuad2DFiniteElement(), H1_TetrahedronElement::H1_TetrahedronElement(), H1_TriangleElement::H1_TriangleElement(), L2_TetrahedronElement::L2_TetrahedronElement(), L2_TriangleElement::L2_TriangleElement(), ND_TetrahedronElement::ND_TetrahedronElement(), ND_TriangleElement::ND_TriangleElement(), RT_TetrahedronElement::RT_TetrahedronElement(), RT_TriangleElement::RT_TriangleElement(), and TestInversion().
| void DenseMatrix::Lump | ( | ) |
Definition at line 1494 of file densemat.cpp.
References Height(), and Width().
Referenced by LumpedIntegrator::AssembleElementMatrix().
Matrix vector multiplication.
Implements Operator.
Definition at line 135 of file densemat.cpp.
References height, mfem_error(), Mult, Vector::Size(), and Operator::size.
| void DenseMatrix::Mult | ( | const double * | x, |
| double * | y | ||
| ) | const |
Matrix vector multiplication.
Definition at line 124 of file densemat.cpp.
References height, and Operator::size.
Multiply a vector with the transpose matrix.
Reimplemented from Operator.
Definition at line 160 of file densemat.cpp.
References height, mfem_error(), Operator::size, and Vector::Size().
Referenced by DiffusionIntegrator::ComputeElementFlux(), GridFunction::ComputeH1Error(), GridFunction::ComputeW11Error(), GridFunction::GetDerivative(), Mesh::GetElementSize(), GridFunction::GetVectorValue(), and GridFunction::GetVectorValues().
| void DenseMatrix::Neg | ( | ) |
(*this) = -(*this)
Definition at line 333 of file densemat.cpp.
References data, Height(), and Operator::Size().
| void DenseMatrix::Norm2 | ( | double * | v | ) | const |
Take the 2-norm of the columns of A and store in v.
Definition at line 437 of file densemat.cpp.
References Height(), and Operator::Size().
Referenced by GridFunction::ComputeL1Error(), GridFunction::ComputeL2Error(), and GridFunction::ComputeMaxError().
| double & DenseMatrix::operator() | ( | int | i, |
| int | j | ||
| ) | [inline] |
Returns reference to a_{ij}. Index i, j = 0 .. size-1.
Definition at line 402 of file densemat.hpp.
References data, height, mfem_error(), and Operator::size.
| const double & DenseMatrix::operator() | ( | int | i, |
| int | j | ||
| ) | const [inline] |
Returns constant reference to a_{ij}. Index i, j = 0 .. size-1.
Definition at line 412 of file densemat.hpp.
References data, height, mfem_error(), and Operator::size.
| double DenseMatrix::operator* | ( | const DenseMatrix & | m | ) | const |
Matrix inner product: tr(A^t B)
Definition at line 145 of file densemat.cpp.
References data, height, Height(), mfem_error(), Operator::size, and Width().
| DenseMatrix & DenseMatrix::operator*= | ( | double | c | ) |
Definition at line 324 of file densemat.cpp.
References data, Height(), and Operator::Size().
| DenseMatrix & DenseMatrix::operator+= | ( | DenseMatrix & | m | ) |
Definition at line 302 of file densemat.cpp.
References height, and Operator::size.
| DenseMatrix & DenseMatrix::operator-= | ( | DenseMatrix & | m | ) |
Definition at line 313 of file densemat.cpp.
References height, and Operator::size.
| DenseMatrix & DenseMatrix::operator= | ( | double | c | ) |
Sets the matrix elements equal to constant c.
Definition at line 280 of file densemat.cpp.
References data, Height(), and Operator::Size().
| DenseMatrix & DenseMatrix::operator= | ( | const DenseMatrix & | m | ) |
Sets the matrix size and elements equal to those of m.
Definition at line 289 of file densemat.cpp.
References data, height, SetSize(), and Operator::size.
| void DenseMatrix::Print | ( | ostream & | out = cout, |
| int | width = 4 |
||
| ) | const [virtual] |
Prints matrix to stream out.
Reimplemented from Matrix.
Definition at line 1733 of file densemat.cpp.
References height, and Operator::size.
| void DenseMatrix::PrintT | ( | ostream & | out = cout, |
| int | width = 4 |
||
| ) | const [virtual] |
Prints the transpose matrix to stream out.
Definition at line 1751 of file densemat.cpp.
References height, and Operator::size.
| void DenseMatrix::SetSize | ( | int | h, |
| int | w | ||
| ) |
If the matrix is not a matrix of size (h x w) then recreate it.
Definition at line 95 of file densemat.cpp.
References data, height, Height(), Operator::size, and Operator::Size().
| void DenseMatrix::SetSize | ( | int | s | ) |
If the matrix is not a square matrix of size s then recreate it.
Definition at line 77 of file densemat.cpp.
References data, height, Height(), Operator::size, and Operator::Size().
Referenced by ElasticityIntegrator::AssembleElementMatrix(), VectorDiffusionIntegrator::AssembleElementMatrix(), DivDivIntegrator::AssembleElementMatrix(), VectorFEMassIntegrator::AssembleElementMatrix(), CurlCurlIntegrator::AssembleElementMatrix(), VectorMassIntegrator::AssembleElementMatrix(), ConvectionIntegrator::AssembleElementMatrix(), MassIntegrator::AssembleElementMatrix(), DiffusionIntegrator::AssembleElementMatrix(), VectorDivergenceIntegrator::AssembleElementMatrix2(), VectorFEMassIntegrator::AssembleElementMatrix2(), DerivativeIntegrator::AssembleElementMatrix2(), VectorFECurlIntegrator::AssembleElementMatrix2(), VectorFEDivergenceIntegrator::AssembleElementMatrix2(), MassIntegrator::AssembleElementMatrix2(), VectorFEDomainLFIntegrator::AssembleRHSElementVect(), VectorFiniteElement::CalcVShape_ND(), Mesh::CheckDisplacements(), DiffusionIntegrator::ComputeElementFlux(), BilinearForm::ComputeElementMatrix(), DiffusionIntegrator::ComputeFluxEnergy(), GridFunction::ComputeH1Error(), GridFunction::ComputeW11Error(), CopyCols(), CopyMN(), CopyRows(), DenseMatrixEigensystem::DenseMatrixEigensystem(), Diag(), dsyev_Eigensystem(), dsyevr_Eigensystem(), MatrixFunctionCoefficient::Eval(), VectorCoefficient::Eval(), NURBSPatch::Get3DRotationMatrix(), Mesh::GetBdrElementTransformation(), Mesh::GetBdrPointMatrix(), GridFunction::GetDerivative(), Mesh::GetElementTransformation(), Mesh::GetFaceTransformation(), Mesh::GetFineElemTrans(), GridFunction::GetGradients(), Mesh::GetLocalQuadToHexTransformation(), Mesh::GetLocalSegToQuadTransformation(), Mesh::GetLocalSegToTriTransformation(), Mesh::GetLocalTriToTetTransformation(), Geometry::GetPerfPointMat(), Mesh::GetPointMatrix(), GridFunction::GetVectorFieldValues(), GridFunction::GetVectorGradient(), GridFunction::GetVectorGradientHat(), GridFunction::GetVectorValues(), H1_TetrahedronElement::H1_TetrahedronElement(), H1_TriangleElement::H1_TriangleElement(), FiniteElementSpace::H2L_GlobalRestrictionMatrix(), IsoparametricTransformation::Jacobian(), L2_TetrahedronElement::L2_TetrahedronElement(), L2_TriangleElement::L2_TriangleElement(), LagrangeHexFiniteElement::LagrangeHexFiniteElement(), ND_TetrahedronElement::ND_TetrahedronElement(), ND_TriangleElement::ND_TriangleElement(), operator=(), Mesh::PrintCharacteristics(), NodalFiniteElement::Project(), VectorFiniteElement::Project_ND(), VectorFiniteElement::Project_RT(), VectorFiniteElement::ProjectCurl_RT(), NodalFiniteElement::ProjectDiv(), NodalFiniteElement::ProjectGrad(), VectorFiniteElement::ProjectGrad_ND(), VectorFiniteElement::ProjectGrad_RT(), RT_TetrahedronElement::RT_TetrahedronElement(), RT_TriangleElement::RT_TriangleElement(), IsoparametricTransformation::Transform(), and Transpose().
| void DenseMatrix::SingularValues | ( | Vector & | sv | ) | const |
Definition at line 699 of file densemat.cpp.
References data, Height(), mfem_error(), Vector::SetSize(), and Width().
| void DenseMatrix::Symmetrize | ( | ) |
(*this) = 1/2 ((*this) + (*this)^t)
Definition at line 1479 of file densemat.cpp.
References Height(), mfem_error(), and Width().
| void DenseMatrix::TestInversion | ( | ) |
Invert and print the numerical conditioning of the inversion.
Definition at line 1769 of file densemat.cpp.
References FNorm(), Invert(), Mult, and Operator::size.
| double DenseMatrix::Trace | ( | ) | const |
Trace of a square matrix.
Definition at line 203 of file densemat.cpp.
References Height(), mfem_error(), Operator::size, and Width().
| void DenseMatrix::Transpose | ( | DenseMatrix & | A | ) |
(*this) = A^t
Definition at line 1470 of file densemat.cpp.
References Height(), SetSize(), and Operator::Size().
| void DenseMatrix::Transpose | ( | ) |
(*this) = (*this)^t
Definition at line 1448 of file densemat.cpp.
References Height(), and Operator::Size().
Referenced by TransposeIntegrator::AssembleElementMatrix(), and TransposeIntegrator::AssembleElementMatrix2().
| void DenseMatrix::UseExternalData | ( | double * | d, |
| int | h, | ||
| int | w | ||
| ) | [inline] |
Definition at line 49 of file densemat.hpp.
References data, height, and Operator::size.
Referenced by DenseTensor::SetSize().
| double DenseMatrix::Weight | ( | ) | const |
Definition at line 250 of file densemat.cpp.
References data, Det(), Height(), mfem_error(), and Width().
Referenced by IsoparametricTransformation::Weight().
| int DenseMatrix::Width | ( | ) | const [inline] |
Returns the width of the matrix.
Definition at line 64 of file densemat.hpp.
References Operator::size.
Referenced by AddMatrix(), AddMult_a_AAt(), AddMult_a_VVt(), AddMultABt(), AddMultADAt(), AddMultVWt(), AdjustDofDirection(), CalcEigenvalues(), CalcInverse(), CalcSingularvalue(), VectorFiniteElement::CalcVShape_ND(), GridFunction::ComputeL1Error(), GridFunction::ComputeL2Error(), GridFunction::ComputeMaxError(), DenseMatrixSVD::DenseMatrixSVD(), Det(), dsyev_Eigensystem(), dsyevr_Eigensystem(), SparseMatrix::EliminateRowColMultipleRHS(), DenseMatrixSVD::Eval(), DenseTensor::GetData(), GradToCurl(), GradToDiv(), IsoparametricTransformation::Jacobian(), Lump(), Mult_a_AAt(), MultABt(), MultADAt(), MultAtB(), operator*(), Mesh::PrintVTK(), GridFunction::SaveVTK(), SingularValues(), DenseTensor::SizeJ(), Symmetrize(), Trace(), and Weight().
friend class DenseMatrixInverse [friend] |
Definition at line 24 of file densemat.hpp.
Referenced by Inverse().
friend class DenseTensor [friend] |
Definition at line 18 of file densemat.hpp.
| void Mult | ( | const DenseMatrix & | b, |
| const DenseMatrix & | c, | ||
| DenseMatrix & | a | ||
| ) | [friend] |
Matrix matrix multiplication. A = B * C.
Definition at line 1804 of file densemat.cpp.
Referenced by ConvectionIntegrator::AssembleElementMatrix(), RT_TetrahedronElement::CalcDivShape(), RT_TriangleElement::CalcDivShape(), L2_TetrahedronElement::CalcShape(), L2_TriangleElement::CalcShape(), H1_TetrahedronElement::CalcShape(), H1_TriangleElement::CalcShape(), GaussQuad2DFiniteElement::CalcShape(), DiffusionIntegrator::ComputeElementFlux(), DiffusionIntegrator::ComputeFluxEnergy(), GridFunction::GetElementAverages(), VectorFiniteElement::LocalInterpolation_ND(), VectorFiniteElement::LocalInterpolation_RT(), Mult(), VectorFiniteElement::Project_ND(), Revolve3D(), NURBSPatch::Rotate3D(), TestInversion(), and IsoparametricTransformation::Transform().
double* DenseMatrix::data [private] |
Definition at line 22 of file densemat.hpp.
Referenced by AddMatrix(), AddToVector(), CalcEigenvalues(), CalcSingularvalue(), ClearExternalData(), CopyMN(), CopyMNt(), Data(), DenseMatrix(), Det(), Diag(), FNorm(), GetColumn(), GetColumnReference(), GetFromVector(), GradToDiv(), Invert(), Mult(), Neg(), operator()(), DenseTensor::operator()(), operator*(), operator*=(), operator=(), SetSize(), SingularValues(), UseExternalData(), Weight(), and ~DenseMatrix().
int DenseMatrix::height [private] |
Definition at line 21 of file densemat.hpp.
Referenced by AddMult(), AddToVector(), ClearExternalData(), DenseMatrix(), DenseMatrixInverse::Factor(), GetColumnReference(), GetFromVector(), GradToDiv(), Height(), InnerProduct(), Mult(), Mult(), MultTranspose(), operator()(), operator*(), operator+=(), operator-=(), operator=(), Print(), PrintT(), SetSize(), and UseExternalData().
1.7.4