28 using AddIntegratorMarkersFn =
33 std::map<BilinearFormIntegrator*, const IntegrationRule*> ir_map;
41 GetIntegratorsFn get_integrators,
42 AddIntegratorFn add_integrator,
50 GetIntegratorsFn get_integrators,
51 GetMarkersFn get_markers,
52 AddIntegratorMarkersFn add_integrator_marker,
53 AddIntegratorFn add_integrator,
58 void ResetIntegrationRules(GetIntegratorsFn get_integrators);
60 static inline int absdof(
int i) {
return i < 0 ? -1-i : i; }
202template <
typename SolverType>
232 template <
typename...
Args>
244 template <
typename...
Args>
246 :
LORSolver(*lor_.GetAssembledSystem(), lor_, args...) { }
@ GaussLobatto
Closed type.
Efficient batched assembly of LOR discretizations on device.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
The Auxiliary-space Divergence Solver in hypre.
The Auxiliary-space Maxwell Solver in hypre.
Wrapper for hypre's ParCSR matrix class.
Class for an integration rule - an Array of IntegrationPoint.
Container class for integration rules.
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
LORBase(FiniteElementSpace &fes_ho_, int ref_type_)
Construct the LORBase object for the given FE space and refinement type.
int GetLOROrder() const
Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
void ConstructLocalDofPermutation(Array< int > &perm_) const
bool HasSameDofNumbering() const
void ConstructDofPermutation() const
void AssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho.
FESpaceType GetFESpaceType() const
Returns the type of finite element space: H1, ND, RT or L2.
FiniteElementCollection * fec
void SetupProlongationAndRestriction()
class BatchedLORAssembly * batched_lor
FiniteElementSpace & fes_ho
virtual void FormLORSpace()=0
Construct the LOR space (overridden for serial and parallel versions).
const Array< int > & GetDofPermutation() const
Returns the permutation that maps LOR DOFs to high-order DOFs.
FiniteElementSpace & GetFESpace() const
Returns the low-order refined finite element space.
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system (const version).
void LegacyAssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho using the legacy method.
Create and assemble a low-order refined version of a BilinearForm.
LORDiscretization(BilinearForm &a_ho, const Array< int > &ess_tdof_list, int ref_type=BasisType::GaussLobatto)
Construct the low-order refined version of a_ho using the given list of essential DOFs.
void FormLORSpace() override
Construct the LOR space (overridden for serial and parallel versions).
SparseMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a SparseMatrix.
OperatorHandle A
The assembled system matrix.
OperatorHandle A
The assembled system matrix.
Represents a solver of type SolverType created using the low-order refined version of the given Bilin...
LORSolver(LORBase &lor_, Args &&... args)
Create a solver of type SolverType using the assembled LOR operator represented by lor_.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
LORSolver(const Operator &op, LORBase &lor_, Args &&... args)
Create a solver of type SolverType using Operator op and arguments args.
const SolverType & GetSolver() const
Access the underlying solver.
SolverType & GetSolver()
Access the underlying solver.
const LORBase & GetLOR() const
Access the LOR discretization object.
LORSolver(ParBilinearForm &a_ho, const Array< int > &ess_tdof_list, int ref_type=BasisType::GaussLobatto)
Create a solver of type SolverType, formed using the assembled HypreParMatrix of the LOR version of a...
LORSolver(BilinearForm &a_ho, const Array< int > &ess_tdof_list, int ref_type=BasisType::GaussLobatto)
Create a solver of type SolverType, formed using the assembled SparseMatrix of the LOR version of a_h...
Pointer to an Operator of a specified type.
int width
Dimension of the input / number of columns in the matrix.
int height
Dimension of the output / number of rows in the matrix.
Abstract parallel finite element space.
Create and assemble a low-order refined version of a ParBilinearForm.
ParFiniteElementSpace & GetParFESpace() const
Return the LOR ParFiniteElementSpace.
ParLORDiscretization(ParBilinearForm &a_ho, const Array< int > &ess_tdof_list, int ref_type=BasisType::GaussLobatto)
Construct the low-order refined version of a_ho using the given list of essential DOFs.
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
void FormLORSpace() override
Construct the LOR space (overridden for serial and parallel versions).