MFEM v4.7.0
Finite element discretization library
|
#include <elasticity_operator.hpp>
Public Member Functions | |
ElasticityOperator (ParMesh &mesh, const int order) | |
virtual void | Mult (const Vector &U, Vector &Y) const override |
Compute the residual Y = R(U) representing the elasticity equation with a material model chosen by calling SetMaterial. | |
Operator & | GetGradient (const Vector &U) const override |
Get the Gradient object. | |
void | GradientMult (const Vector &dX, Vector &Y) const |
Multiply the linearization of the residual R(U) wrt to the current state U by a perturbation dX. | |
void | AssembleGradientDiagonal (Vector &Ke_diag, Vector &K_diag_local, Vector &K_diag) const |
Assemble the linearization of the residual K = dR(U)/dU. | |
~ElasticityOperator () | |
template<typename material_type > | |
void | SetMaterial (const material_type &material) |
Set the material type. | |
void | SetEssentialAttributes (const Array< int > attr) |
Set the essential attributes which mark degrees of freedom for the solving process. | |
void | SetPrescribedDisplacement (const Array< int > attr) |
Set the attributes which mark the degrees of freedom that have a fixed displacement. | |
const Array< int > & | GetPrescribedDisplacementTDofs () |
Return the T-vector degrees of freedom that have been marked as displaced. | |
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. | |
Operator (int s=0) | |
Construct a square Operator with given size s (default 0). | |
Operator (int h, int w) | |
Construct an Operator with the given height (output size) and width (input size). | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows(). | |
int | NumRows () const |
Get the number of rows (size of output) of the Operator. Synonym with Height(). | |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols(). | |
int | NumCols () const |
Get the number of columns (size of input) of the Operator. Synonym with Width(). | |
virtual MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. | |
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. | |
virtual void | AddMult (const Vector &x, Vector &y, const real_t a=1.0) const |
Operator application: y+=A(x) (default) or y+=a*A(x) . | |
virtual void | AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const |
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x) . | |
virtual void | ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const |
Operator application on a matrix: Y=A(X) . | |
virtual void | ArrayMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y) const |
Action of the transpose operator on a matrix: Y=A^t(X) . | |
virtual void | ArrayAddMult (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const |
Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X) . | |
virtual void | ArrayAddMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const |
Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X) . | |
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. | |
virtual const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. | |
virtual const Operator * | GetOutputProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. | |
virtual const Operator * | GetOutputRestrictionTranspose () const |
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. | |
virtual const Operator * | GetOutputRestriction () const |
Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. | |
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. | |
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. | |
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(). | |
void | FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this square operator. | |
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). | |
void | FormDiscreteOperator (Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator. | |
void | PrintMatlab (std::ostream &out, int n, int m=0) const |
Prints operator with input size n and output size m in Matlab format. | |
virtual void | PrintMatlab (std::ostream &out) const |
Prints operator in Matlab format. | |
virtual | ~Operator () |
Virtual destructor. | |
Type | GetType () const |
Return the type ID of the Operator class. | |
Public Attributes | |
ParMesh & | mesh_ |
const int | order_ |
Polynomial order of the FE space. | |
const int | dim_ |
const int | vdim_ |
const int | ne_ |
Number of elements in the mesh (rank local) | |
H1_FECollection | h1_fec_ |
H1 finite element collection. | |
ParFiniteElementSpace | h1_fes_ |
H1 finite element space. | |
IntegrationRule * | ir_ = nullptr |
int | d1d_ |
Number of degrees of freedom in 1D. | |
int | q1d_ |
Number of quadrature points in 1D. | |
const Operator * | h1_element_restriction_ |
const Operator * | h1_prolongation_ |
Array< int > | ess_tdof_list_ |
Array< int > | displaced_tdof_list_ |
Operator * | gradient_ |
const GeometricFactors * | geometric_factors_ |
const DofToQuad * | maps_ |
Vector | X_local_ |
Input state L-vector. | |
Vector | X_el_ |
Input state E-vector. | |
Vector | Y_local_ |
Output state L-vector. | |
Vector | Y_el_ |
Output state E-Vector. | |
Vector | current_state |
Vector | cstate_local_ |
Vector | cstate_el_ |
Vector | dX_ess_ |
bool | use_cache_ = true |
bool | recompute_cache_ = false |
Vector | dsigma_cache_ |
std::function< void(const int, const Array< real_t > &, const Array< real_t > &, const Array< real_t > &, const Vector &, const Vector &, const Vector &, Vector &)> | element_apply_kernel_wrapper |
Wrapper for the application of the residual R(U). | |
std::function< void(const int, const Array< real_t > &, const Array< real_t > &, const Array< real_t > &, const Vector &, const Vector &, const Vector &, Vector &, const Vector &)> | element_apply_gradient_kernel_wrapper |
Wrapper for the application of the gradient of the residual. | |
std::function< void(const int, const Array< real_t > &, const Array< real_t > &, const Array< real_t > &, const Vector &, const Vector &, const Vector &, Vector &)> | element_kernel_assemble_diagonal_wrapper |
Wrapper for the assembly of the gradient on each diagonal element. | |
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() | |
void | FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout) |
see FormRectangularSystemOperator() | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". | |
Protected Attributes inherited from mfem::Operator | |
int | height |
Dimension of the output / number of rows in the matrix. | |
int | width |
Dimension of the input / number of columns in the matrix. | |
Definition at line 20 of file elasticity_operator.hpp.
mfem::ElasticityOperator::ElasticityOperator | ( | ParMesh & | mesh, |
const int | order ) |
Definition at line 18 of file elasticity_operator.cpp.
mfem::ElasticityOperator::~ElasticityOperator | ( | ) |
Definition at line 211 of file elasticity_operator.cpp.
void mfem::ElasticityOperator::AssembleGradientDiagonal | ( | Vector & | Ke_diag, |
Vector & | K_diag_local, | ||
Vector & | K_diag ) const |
Assemble the linearization of the residual K = dR(U)/dU.
This method needs three input vectors which also act as output vectors. They don't have to be the right size on the first call, but it is advised that memory is kept alive during successive call. The data layout of the outputs will be
Ke_diag: dofs x dofs x dofs x dim x ne x dim
K_diag_local: width(H1_Restriction) x dim
K_diag: width(H1_Prolongation) x dim
This data layout is needed due to the Ordering::byNODES. See method implementation comments for more details. The output K_diag has modified entries when essential boundaries are defined. Each essential dof row and column are set to zero with it's diagonal entry set to 1.
Ke_diag | |
K_diag_local | |
K_diag |
Definition at line 150 of file elasticity_operator.cpp.
Get the Gradient object.
Update and cache the state vector U, used to compute the linearization dR(U)/dU.
U |
Reimplemented from mfem::Operator.
Definition at line 100 of file elasticity_operator.cpp.
|
inline |
Return the T-vector degrees of freedom that have been marked as displaced.
Definition at line 300 of file elasticity_operator.hpp.
Multiply the linearization of the residual R(U) wrt to the current state U by a perturbation dX.
Y = dR(U)/dU * dX = K(U) dX
dX | |
Y |
Definition at line 110 of file elasticity_operator.cpp.
Compute the residual Y = R(U) representing the elasticity equation with a material model chosen by calling SetMaterial.
The output vector Y has essential degrees of freedom applied by setting them to zero. This ensures R(U)_i = 0 being satisfied for each essential dof i.
U | U |
Y | Residual R(U) |
Implements mfem::Operator.
Definition at line 74 of file elasticity_operator.cpp.
|
inline |
Set the essential attributes which mark degrees of freedom for the solving process.
Can be either a fixed boundary or a prescribed displacement.
attr |
Definition at line 278 of file elasticity_operator.hpp.
|
inline |
Set the material type.
This method sets the material type by instantiating the kernels with a material_type object.
material_type |
[in] | material |
Definition at line 179 of file elasticity_operator.hpp.
|
inline |
Set the attributes which mark the degrees of freedom that have a fixed displacement.
[in] | attr |
Definition at line 289 of file elasticity_operator.hpp.
|
mutable |
Definition at line 124 of file elasticity_operator.hpp.
|
mutable |
Definition at line 123 of file elasticity_operator.hpp.
|
mutable |
Cached current state. Used to determine the state on which to compute the linearization on during the Newton method.
Definition at line 122 of file elasticity_operator.hpp.
int mfem::ElasticityOperator::d1d_ |
Number of degrees of freedom in 1D.
Definition at line 102 of file elasticity_operator.hpp.
const int mfem::ElasticityOperator::dim_ |
Definition at line 91 of file elasticity_operator.hpp.
Array<int> mfem::ElasticityOperator::displaced_tdof_list_ |
Definition at line 108 of file elasticity_operator.hpp.
Vector mfem::ElasticityOperator::dsigma_cache_ |
Definition at line 134 of file elasticity_operator.hpp.
|
mutable |
Temporary vector for the perturbation of the solution with essential boundaries eliminated. Defined as a T-vector.
Definition at line 127 of file elasticity_operator.hpp.
std::function<void(const int, const Array<real_t> &, const Array<real_t> &, const Array<real_t> &, const Vector &, const Vector &, const Vector &, Vector &, const Vector &)> mfem::ElasticityOperator::element_apply_gradient_kernel_wrapper |
Wrapper for the application of the gradient of the residual.
K(U) dX = dR(U)/dU dX
Definition at line 157 of file elasticity_operator.hpp.
std::function<void(const int, const Array<real_t> &, const Array<real_t> &, const Array<real_t> &, const Vector &, const Vector &, const Vector &, Vector &)> mfem::ElasticityOperator::element_apply_kernel_wrapper |
Wrapper for the application of the residual R(U).
The wrapper is used in SetMaterial to instantiate the chosen kernel and erase the material type kernel. This is purely an interface design choice and could be replaced by an abstract base class for the material including virtual function calls.
Definition at line 147 of file elasticity_operator.hpp.
std::function<void(const int, const Array<real_t> &, const Array<real_t> &, const Array<real_t> &, const Vector &, const Vector &, const Vector &, Vector &)> mfem::ElasticityOperator::element_kernel_assemble_diagonal_wrapper |
Wrapper for the assembly of the gradient on each diagonal element.
Ke_ii(U) = dRe_ii(U)/dU
Definition at line 167 of file elasticity_operator.hpp.
Array<int> mfem::ElasticityOperator::ess_tdof_list_ |
Definition at line 107 of file elasticity_operator.hpp.
const GeometricFactors* mfem::ElasticityOperator::geometric_factors_ |
Definition at line 110 of file elasticity_operator.hpp.
Operator* mfem::ElasticityOperator::gradient_ |
Definition at line 109 of file elasticity_operator.hpp.
const Operator* mfem::ElasticityOperator::h1_element_restriction_ |
Definition at line 105 of file elasticity_operator.hpp.
H1_FECollection mfem::ElasticityOperator::h1_fec_ |
H1 finite element collection.
Definition at line 96 of file elasticity_operator.hpp.
ParFiniteElementSpace mfem::ElasticityOperator::h1_fes_ |
H1 finite element space.
Definition at line 98 of file elasticity_operator.hpp.
const Operator* mfem::ElasticityOperator::h1_prolongation_ |
Definition at line 106 of file elasticity_operator.hpp.
IntegrationRule* mfem::ElasticityOperator::ir_ = nullptr |
Definition at line 100 of file elasticity_operator.hpp.
const DofToQuad* mfem::ElasticityOperator::maps_ |
Definition at line 111 of file elasticity_operator.hpp.
ParMesh& mfem::ElasticityOperator::mesh_ |
Definition at line 88 of file elasticity_operator.hpp.
const int mfem::ElasticityOperator::ne_ |
Number of elements in the mesh (rank local)
Definition at line 94 of file elasticity_operator.hpp.
const int mfem::ElasticityOperator::order_ |
Polynomial order of the FE space.
Definition at line 90 of file elasticity_operator.hpp.
int mfem::ElasticityOperator::q1d_ |
Number of quadrature points in 1D.
Definition at line 104 of file elasticity_operator.hpp.
|
mutable |
Definition at line 133 of file elasticity_operator.hpp.
bool mfem::ElasticityOperator::use_cache_ = true |
Flag to enable caching of the gradient. If enabled, during linear iterations the operator only applies the gradient on each quadrature point rather than recompute the action.
Definition at line 132 of file elasticity_operator.hpp.
const int mfem::ElasticityOperator::vdim_ |
Definition at line 92 of file elasticity_operator.hpp.
|
mutable |
Input state E-vector.
Definition at line 115 of file elasticity_operator.hpp.
|
mutable |
Input state L-vector.
Definition at line 113 of file elasticity_operator.hpp.
|
mutable |
Output state E-Vector.
Definition at line 119 of file elasticity_operator.hpp.
|
mutable |
Output state L-vector.
Definition at line 117 of file elasticity_operator.hpp.