MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t > Class Template Reference

Templated bilinear form class, cf. bilinearform.?pp. More...

#include <tbilinearform.hpp>

Inheritance diagram for mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >:
[legend]
Collaboration diagram for mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >:
[legend]

Classes

struct  S_spec
 Contains matrix sizes, type of kernel (ElementMatrix is templated on a kernel, e.g. ElementMatrix::Compute may be AssembleGradGrad()). More...
 
struct  T_result
 

Public Types

typedef impl_traits_t impl_traits_type
 
- 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...
 

Public Member Functions

 TBilinearForm (const IntegratorType &integ, const FiniteElementSpace &sol_fes)
 
virtual ~TBilinearForm ()
 
virtual const OperatorGetProlongation () const
 Get the input finite element space prolongation matrix.
 
virtual const OperatorGetRestriction () const
 Get the input finite element space restriction matrix.
 
virtual void Mult (const Vector &x, Vector &y) const
 Operator application: y=A(x).
 
void MultUnassembled (const Vector &x, Vector &y) const
 
void Assemble ()
 Partial assembly of quadrature point data.
 
MFEM_ALWAYS_INLINE void ElementAddMultAssembled (int el, solFieldEval &solFEval) const
 
void MultAssembled (const Vector &x, Vector &y) const
 
void TestElementwiseExtractAssemble (const Vector &x, Vector &y) const
 
void SerializeNodes (Vector &sNodes) const
 
void AssembleFromSerializedNodes (const Vector &sNodes)
 Partial assembly from "serialized" nodes.
 
void Serialize (const Vector &x, Vector &sx) const
 
void MultAssembledSerialized (const Vector &sx, Vector &sy) const
 serialized vector sx --> serialized vector 'sy'
 
void AssembleMatrix (SparseMatrix &M) const
 Assemble the operator in a SparseMatrix.
 
void AssembleMatrix (DenseTensor &M) const
 Assemble element matrices and store them as a DenseTensor object.
 
void AssembleBilinearForm (BilinearForm &a) const
 Assemble element matrices and add them to the bilinear form.
 
void AddMult (DenseTensor &M, const Vector &x, Vector &y) const
 Multiplication using assembled element matrices stored as a DenseTensor.
 
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).
 
- 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 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 OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error.
 
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 OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity.
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators.
 
virtual const OperatorGetOutputRestriction () 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.
 

Protected Types

typedef complex_t complex_type
 
typedef real_t real_type
 
typedef meshType::FE_type meshFE_type
 
typedef ShapeEvaluator< meshFE_type, IR, real_tmeshShapeEval
 
typedef solFESpace::FE_type solFE_type
 
typedef ShapeEvaluator< solFE_type, IR, real_tsolShapeEval
 
typedef solVecLayout_t solVecLayout_type
 
typedef impl_traits_t::vcomplex_t vcomplex_t
 
typedef impl_traits_t::vreal_t vreal_t
 
typedef kernel_t::template CoefficientEval< IR, coeff_t, impl_traits_t >::Type coeff_eval_t
 
typedef TElementTransformation< meshType, IR, real_tTrans_t
 
typedef FieldEvaluator< solFESpace, solVecLayout_t, IR, complex_t, real_tsolFieldEval
 
IntegratorType defines several internal types
typedef IntegratorType integ_t
 
typedef integ_t::coefficient_type coeff_t
 coeff_t might be TConstantCoefficient or TFunctionCoefficient, for example
 
typedef integ_t::template kernel< sdim, dim, vcomplex_t >::type kernel_t
 kernel_t may be TDiffusionKernel or TMassKernel
 
typedef kernel_t::template p_asm_data< qpts >::type p_assembled_t
 p_assembled_t is something like a TTensor or TMatrix for partial assembly
 
typedef kernel_t::template f_asm_data< qpts >::type f_assembled_t
 f_assembled_t is something like a TTensor or TMatrix for full assembly
 

Protected Attributes

meshType mesh
 
meshShapeEval meshEval
 
solFE_type sol_fe
 
solShapeEval solEval
 
solFESpace solFES
 
solVecLayout_t solVecLayout
 
IR int_rule
 
coeff_t coeff
 
Memory< p_assembled_tassembled_data
 
const FiniteElementSpacein_fes
 
- 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.
 

Static Protected Attributes

static const int dim = meshType::dim
 
static const int sdim = meshType::space_dim
 
static const int dofs = solFE_type::dofs
 
static const int vdim = solVecLayout_t::vec_dim
 
static const int qpts = IR::qpts
 
static const int AB = impl_traits_t::align_bytes
 
static const int SS = impl_traits_t::simd_size
 
static const int BE = impl_traits_t::batch_size
 
static const int TE = SS*BE
 

Additional Inherited Members

- 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()
 
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".
 

Detailed Description

template<typename meshType, typename solFESpace, typename IR, typename IntegratorType, typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
class mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >

Templated bilinear form class, cf. bilinearform.?pp.

Template Parameters
meshTypetypically TMesh, which is templated on FE type
solFESpaceeg. H1_FiniteElementSpace
IRintegration rule, typically TIntegrationRule, which is further templated on element geometry
IntegratorTypetypically a TIntegrator, which is templated on a kernel, eg. TDiffusionKernel or TMassKernel. This describes what actual problem you solve.
solVecLayout_tdescribes how degrees of freedom are laid out, scalar or vector, column/row major, etc.
complex_tdata type for solution dofs
real_tdata type for mesh nodes, solution basis, and mesh basis

Definition at line 46 of file tbilinearform.hpp.

Member Typedef Documentation

◆ coeff_eval_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef kernel_t::template CoefficientEval<IR,coeff_t,impl_traits_t>::Type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::coeff_eval_t
protected

Definition at line 88 of file tbilinearform.hpp.

◆ coeff_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef integ_t::coefficient_type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::coeff_t
protected

coeff_t might be TConstantCoefficient or TFunctionCoefficient, for example

Definition at line 78 of file tbilinearform.hpp.

◆ complex_type

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef complex_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::complex_type
protected

Definition at line 52 of file tbilinearform.hpp.

◆ f_assembled_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef kernel_t::template f_asm_data<qpts>::type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::f_assembled_t
protected

f_assembled_t is something like a TTensor or TMatrix for full assembly

Definition at line 84 of file tbilinearform.hpp.

◆ impl_traits_type

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef impl_traits_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::impl_traits_type

Definition at line 49 of file tbilinearform.hpp.

◆ integ_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef IntegratorType mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::integ_t
protected

Definition at line 76 of file tbilinearform.hpp.

◆ kernel_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef integ_t::template kernel<sdim,dim,vcomplex_t>::type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::kernel_t
protected

kernel_t may be TDiffusionKernel or TMassKernel

Definition at line 80 of file tbilinearform.hpp.

◆ meshFE_type

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef meshType::FE_type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::meshFE_type
protected

Definition at line 55 of file tbilinearform.hpp.

◆ meshShapeEval

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef ShapeEvaluator<meshFE_type,IR,real_t> mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::meshShapeEval
protected

Definition at line 56 of file tbilinearform.hpp.

◆ p_assembled_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef kernel_t::template p_asm_data<qpts>::type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::p_assembled_t
protected

p_assembled_t is something like a TTensor or TMatrix for partial assembly

Definition at line 82 of file tbilinearform.hpp.

◆ real_type

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef real_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::real_type
protected

Definition at line 53 of file tbilinearform.hpp.

◆ solFE_type

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef solFESpace::FE_type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::solFE_type
protected

Definition at line 57 of file tbilinearform.hpp.

◆ solFieldEval

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef FieldEvaluator<solFESpace,solVecLayout_t,IR, complex_t,real_t> mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::solFieldEval
protected

Definition at line 100 of file tbilinearform.hpp.

◆ solShapeEval

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef ShapeEvaluator<solFE_type,IR,real_t> mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::solShapeEval
protected

Definition at line 58 of file tbilinearform.hpp.

◆ solVecLayout_type

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef solVecLayout_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::solVecLayout_type
protected

Definition at line 59 of file tbilinearform.hpp.

◆ Trans_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef TElementTransformation<meshType,IR,real_t> mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::Trans_t
protected

Definition at line 91 of file tbilinearform.hpp.

◆ vcomplex_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef impl_traits_t::vcomplex_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::vcomplex_t
protected

Definition at line 71 of file tbilinearform.hpp.

◆ vreal_t

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
typedef impl_traits_t::vreal_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::vreal_t
protected

Definition at line 72 of file tbilinearform.hpp.

Constructor & Destructor Documentation

◆ TBilinearForm()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::TBilinearForm ( const IntegratorType & integ,
const FiniteElementSpace & sol_fes )
inline

Definition at line 130 of file tbilinearform.hpp.

◆ ~TBilinearForm()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
virtual mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::~TBilinearForm ( )
inlinevirtual

Definition at line 148 of file tbilinearform.hpp.

Member Function Documentation

◆ AddMult() [1/2]

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::Operator::AddMult ( const Vector & x,
Vector & y,
const real_t a = 1.0 ) const
virtual

Operator application: y+=A(x) (default) or y+=a*A(x).

Reimplemented from mfem::Operator.

Definition at line 97 of file operator.cpp.

◆ AddMult() [2/2]

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::AddMult ( DenseTensor & M,
const Vector & x,
Vector & y ) const
inline

Multiplication using assembled element matrices stored as a DenseTensor.

Definition at line 598 of file tbilinearform.hpp.

◆ Assemble()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::Assemble ( )
inline

Partial assembly of quadrature point data.

Definition at line 216 of file tbilinearform.hpp.

◆ AssembleBilinearForm()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::AssembleBilinearForm ( BilinearForm & a) const
inline

Assemble element matrices and add them to the bilinear form.

Definition at line 498 of file tbilinearform.hpp.

◆ AssembleFromSerializedNodes()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::AssembleFromSerializedNodes ( const Vector & sNodes)
inline

Partial assembly from "serialized" nodes.

Definition at line 317 of file tbilinearform.hpp.

◆ AssembleMatrix() [1/2]

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::AssembleMatrix ( DenseTensor & M) const
inline

Assemble element matrices and store them as a DenseTensor object.

Definition at line 445 of file tbilinearform.hpp.

◆ AssembleMatrix() [2/2]

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::AssembleMatrix ( SparseMatrix & M) const
inline

Assemble the operator in a SparseMatrix.

Definition at line 397 of file tbilinearform.hpp.

◆ ElementAddMultAssembled()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
MFEM_ALWAYS_INLINE void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::ElementAddMultAssembled ( int el,
solFieldEval & solFEval ) const
inline

Definition at line 243 of file tbilinearform.hpp.

◆ GetProlongation()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
virtual const Operator * mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::GetProlongation ( ) const
inlinevirtual

Get the input finite element space prolongation matrix.

Reimplemented from mfem::Operator.

Definition at line 154 of file tbilinearform.hpp.

◆ GetRestriction()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
virtual const Operator * mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::GetRestriction ( ) const
inlinevirtual

Get the input finite element space restriction matrix.

Reimplemented from mfem::Operator.

Definition at line 157 of file tbilinearform.hpp.

◆ Mult()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
virtual void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::Mult ( const Vector & x,
Vector & y ) const
inlinevirtual

Operator application: y=A(x).

Implements mfem::Operator.

Definition at line 160 of file tbilinearform.hpp.

◆ MultAssembled()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::MultAssembled ( const Vector & x,
Vector & y ) const
inline

Definition at line 257 of file tbilinearform.hpp.

◆ MultAssembledSerialized()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::MultAssembledSerialized ( const Vector & sx,
Vector & sy ) const
inline

serialized vector sx --> serialized vector 'sy'

Definition at line 369 of file tbilinearform.hpp.

◆ MultUnassembled()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::MultUnassembled ( const Vector & x,
Vector & y ) const
inline

Definition at line 173 of file tbilinearform.hpp.

◆ Serialize()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::Serialize ( const Vector & x,
Vector & sx ) const
inline

Definition at line 345 of file tbilinearform.hpp.

◆ SerializeNodes()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::SerializeNodes ( Vector & sNodes) const
inline

Definition at line 293 of file tbilinearform.hpp.

◆ TestElementwiseExtractAssemble()

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
void mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::TestElementwiseExtractAssemble ( const Vector & x,
Vector & y ) const
inline

Definition at line 273 of file tbilinearform.hpp.

Member Data Documentation

◆ AB

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::AB = impl_traits_t::align_bytes
staticprotected

Definition at line 66 of file tbilinearform.hpp.

◆ assembled_data

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
Memory<p_assembled_t> mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::assembled_data
protected

Definition at line 125 of file tbilinearform.hpp.

◆ BE

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::BE = impl_traits_t::batch_size
staticprotected

Definition at line 68 of file tbilinearform.hpp.

◆ coeff

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
coeff_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::coeff
protected

Definition at line 123 of file tbilinearform.hpp.

◆ dim

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::dim = meshType::dim
staticprotected

Definition at line 61 of file tbilinearform.hpp.

◆ dofs

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::dofs = solFE_type::dofs
staticprotected

Definition at line 63 of file tbilinearform.hpp.

◆ in_fes

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const FiniteElementSpace& mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::in_fes
protected

Definition at line 127 of file tbilinearform.hpp.

◆ int_rule

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
IR mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::int_rule
protected

Definition at line 121 of file tbilinearform.hpp.

◆ mesh

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
meshType mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::mesh
protected

Definition at line 113 of file tbilinearform.hpp.

◆ meshEval

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
meshShapeEval mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::meshEval
protected

Definition at line 114 of file tbilinearform.hpp.

◆ qpts

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::qpts = IR::qpts
staticprotected

Definition at line 65 of file tbilinearform.hpp.

◆ sdim

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::sdim = meshType::space_dim
staticprotected

Definition at line 62 of file tbilinearform.hpp.

◆ sol_fe

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
solFE_type mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::sol_fe
protected

Definition at line 116 of file tbilinearform.hpp.

◆ solEval

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
solShapeEval mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::solEval
protected

Definition at line 117 of file tbilinearform.hpp.

◆ solFES

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
solFESpace mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::solFES
mutableprotected

Definition at line 118 of file tbilinearform.hpp.

◆ solVecLayout

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
solVecLayout_t mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::solVecLayout
protected

Definition at line 119 of file tbilinearform.hpp.

◆ SS

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::SS = impl_traits_t::simd_size
staticprotected

Definition at line 67 of file tbilinearform.hpp.

◆ TE

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::TE = SS*BE
staticprotected

Definition at line 69 of file tbilinearform.hpp.

◆ vdim

template<typename meshType , typename solFESpace , typename IR , typename IntegratorType , typename solVecLayout_t = ScalarLayout, typename complex_t = real_t, typename real_t = real_t, typename impl_traits_t = AutoSIMDTraits<complex_t,real_t>>
const int mfem::TBilinearForm< meshType, solFESpace, IR, IntegratorType, solVecLayout_t, complex_t, real_t, impl_traits_t >::vdim = solVecLayout_t::vec_dim
staticprotected

Definition at line 64 of file tbilinearform.hpp.


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