15 #include "../bilinearform.hpp" 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; }
202 template <
typename SolverType>
232 template <
typename...
Args>
244 template <
typename...
Args>
246 :
LORSolver(*lor_.GetAssembledSystem(), lor_, args...) { }
class BatchedLORAssembly * batched_lor
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
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...
Create and assemble a low-order refined version of a ParBilinearForm.
Class for an integration rule - an Array of IntegrationPoint.
The Auxiliary-space Maxwell Solver in hypre.
The Auxiliary-space Divergence Solver in hypre.
OperatorHandle A
The assembled system matrix.
ParFiniteElementSpace & GetParFESpace() const
Return the LOR ParFiniteElementSpace.
void LegacyAssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho using the legacy method.
int GetLOROrder() const
Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
Efficient batched assembly of LOR discretizations on device.
Pointer to an Operator of a specified type.
FiniteElementSpace & fes_ho
Container class for integration rules.
const SolverType & GetSolver() const
Access the underlying solver.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Abstract parallel finite element space.
void FormLORSpace() override
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.
Create and assemble a low-order refined version of a BilinearForm.
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(LORBase &lor_, Args &&... args)
Create a solver of type SolverType using the assembled LOR operator represented by lor_...
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system (const version).
OperatorHandle A
The assembled system matrix.
const LORBase & GetLOR() const
Access the LOR discretization object.
LORBase(FiniteElementSpace &fes_ho_, int ref_type_)
Construct the LORBase object for the given FE space and refinement type.
void SetupProlongationAndRestriction()
virtual void FormLORSpace()=0
Construct the LOR space (overridden for serial and parallel versions).
bool HasSameDofNumbering() const
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.
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
SolverType & GetSolver()
Access the underlying solver.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
void ConstructDofPermutation() const
SparseMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a SparseMatrix.
int height
Dimension of the output / number of rows in the matrix.
FiniteElementCollection * fec
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...
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.
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...
FiniteElementSpace & GetFESpace() const
Returns the low-order refined finite element space.
void FormLORSpace() override
Construct the LOR space (overridden for serial and parallel versions).
void ConstructLocalDofPermutation(Array< int > &perm_) const
Wrapper for hypre's ParCSR matrix class.
int width
Dimension of the input / number of columns in the matrix.
Represents a solver of type SolverType created using the low-order refined version of the given Bilin...