MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lor.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_LOR
13 #define MFEM_LOR
14 
15 #include "bilinearform.hpp"
16 
17 namespace mfem
18 {
19 
20 /// @brief Abstract base class for LORDiscretization and ParLORDiscretization
21 /// classes, which construct low-order refined versions of bilinear forms.
22 class LORBase
23 {
24 private:
25  using GetIntegratorsFn = Array<BilinearFormIntegrator*> *(BilinearForm::*)();
26  using GetMarkersFn = Array<Array<int>*> *(BilinearForm::*)();
27  using AddIntegratorFn = void (BilinearForm::*)(BilinearFormIntegrator*);
28  using AddIntegratorMarkersFn =
30 
31  IntegrationRules irs;
32  const IntegrationRule *ir_el, *ir_face;
33  std::map<BilinearFormIntegrator*, const IntegrationRule*> ir_map;
34 
35  /// Adds all the integrators from the BilinearForm @a a_from to @a a_to. If
36  /// the mesh consists of tensor product elements, temporarily changes the
37  /// integration rules of the integrators to use collocated quadrature for
38  /// better conditioning of the LOR system.
39  void AddIntegrators(BilinearForm &a_from,
40  BilinearForm &a_to,
41  GetIntegratorsFn get_integrators,
42  AddIntegratorFn add_integrator,
43  const IntegrationRule *ir);
44 
45  /// Adds all the integrators from the BilinearForm @a a_from to @a a_to,
46  /// using the same marker lists (e.g. for integrators only applied to
47  /// certain boundary attributes). @sa LORBase::AddIntegrators
48  void AddIntegratorsAndMarkers(BilinearForm &a_from,
49  BilinearForm &a_to,
50  GetIntegratorsFn get_integrators,
51  GetMarkersFn get_markers,
52  AddIntegratorMarkersFn add_integrator_marker,
53  AddIntegratorFn add_integrator,
54  const IntegrationRule *ir);
55 
56  /// Resets the integration rules of the integrators of @a a to their original
57  /// values (after temporarily changing them for LOR assembly).
58  void ResetIntegrationRules(GetIntegratorsFn get_integrators);
59 
60  static inline int absdof(int i) { return i < 0 ? -1-i : i; }
61 
62 protected:
63  enum FESpaceType { H1, ND, RT, L2, INVALID };
64 
71  mutable Array<int> perm;
72 
73  /// Constructs the local DOF (ldof) permutation. In parallel this is used as
74  /// an intermediate step in computing the DOF permutation (see
75  /// ConstructDofPermutation and GetDofPermutation).
76  void ConstructLocalDofPermutation(Array<int> &perm_) const;
77 
78  /// Construct the permutation that maps LOR DOFs to high-order DOFs. See
79  /// GetDofPermutation.
80  void ConstructDofPermutation() const;
81 
82  /// Returns true if the LOR space and HO space have the same DOF numbering
83  /// (H1 or L2 spaces), false otherwise (ND or RT spaces).
84  bool HasSameDofNumbering() const;
85 
86  /// Sets up the prolongation and restriction operators required in the case
87  /// of different DOF numberings (ND or RT spaces) or nonconforming spaces.
89 
90  /// Returns the type of finite element space: H1, ND, RT or L2.
92 
93  /// Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
94  int GetLOROrder() const;
95 
96  /// Assembles the LOR system (used internally by
97  /// LORDiscretization::AssembleSystem and
98  /// ParLORDiscretization::AssembleSystem).
99  void AssembleSystem_(BilinearForm &a_ho, const Array<int> &ess_dofs);
100 
101  LORBase(FiniteElementSpace &fes_ho_);
102 
103 public:
104  /// Returns the assembled LOR system.
105  const OperatorHandle &GetAssembledSystem() const;
106 
107  /// @brief Returns the permutation that maps LOR DOFs to high-order DOFs.
108  ///
109  /// This permutation is constructed the first time it is requested, and then
110  /// is cached. For H1 and L2 finite element spaces (or for nonconforming
111  /// spaces) this is the identity. In these cases, RequiresDofPermutation will
112  /// return false. However, if the DOF permutation is requested, an identity
113  /// permutation will be built and returned.
114  ///
115  /// For vector finite element spaces (ND and RT), the DOF permutation is
116  /// nontrivial. Returns an array @a perm such that, given an index @a i of a
117  /// LOR dof, @a perm[i] is the index of the corresponding HO dof.
118  const Array<int> &GetDofPermutation() const;
119 
120  /// Returns the low-order refined finite element space.
121  FiniteElementSpace &GetFESpace() const { return *fes; }
122 
123  ~LORBase();
124 };
125 
126 /// Create and assemble a low-order refined version of a BilinearForm.
128 {
129 public:
130  /// @brief Construct the low-order refined version of @a a_ho using the given
131  /// list of essential DOFs.
132  ///
133  /// The mesh is refined using the refinement type specified by @a ref_type
134  /// (see Mesh::MakeRefined).
135  LORDiscretization(BilinearForm &a_ho, const Array<int> &ess_tdof_list,
136  int ref_type=BasisType::GaussLobatto);
137 
138  /// @brief Construct a low-order refined version of the FiniteElementSpace @a
139  /// fes_ho.
140  ///
141  /// The mesh is refined using the refinement type specified by @a ref_type
142  /// (see Mesh::MakeRefined).
144  int ref_type=BasisType::GaussLobatto);
145 
146  /// Assembles the LOR system corresponding to @a a_ho.
147  void AssembleSystem(BilinearForm &a_ho, const Array<int> &ess_dofs);
148 
149  /// Return the assembled LOR operator as a SparseMatrix.
151 };
152 
153 #ifdef MFEM_USE_MPI
154 
155 /// Create and assemble a low-order refined version of a ParBilinearForm.
157 {
158 public:
159  /// @brief Construct the low-order refined version of @a a_ho using the given
160  /// list of essential DOFs.
161  ///
162  /// The mesh is refined using the refinement type specified by @a ref_type
163  /// (see ParMesh::MakeRefined).
164  ParLORDiscretization(ParBilinearForm &a_ho, const Array<int> &ess_tdof_list,
165  int ref_type=BasisType::GaussLobatto);
166 
167  /// @brief Construct a low-order refined version of the ParFiniteElementSpace
168  /// @a pfes_ho.
169  ///
170  /// The mesh is refined using the refinement type specified by @a ref_type
171  /// (see ParMesh::MakeRefined).
173  int ref_type=BasisType::GaussLobatto);
174 
175  /// Assembles the LOR system corresponding to @a a_ho.
176  void AssembleSystem(ParBilinearForm &a_ho, const Array<int> &ess_dofs);
177 
178  /// Return the assembled LOR operator as a HypreParMatrix.
180 
181  /// Return the LOR ParFiniteElementSpace.
183 };
184 
185 #endif
186 
187 /// @brief Represents a solver of type @a SolverType created using the low-order
188 /// refined version of the given BilinearForm or ParBilinearForm.
189 ///
190 /// @note To achieve good solver performance, the high-order finite element
191 /// space should use BasisType::GaussLobatto for H1 discretizations, and basis
192 /// pair (BasisType::GaussLobatto, BasisType::IntegratedGLL) for Nedelec and
193 /// Raviart-Thomas elements.
194 template <typename SolverType>
195 class LORSolver : public Solver
196 {
197 protected:
199  bool own_lor = true;
200  SolverType solver;
201  mutable Vector px, py;
202 public:
203  /// @brief Create a solver of type @a SolverType, formed using the assembled
204  /// SparseMatrix of the LOR version of @a a_ho. @see LORDiscretization
205  LORSolver(BilinearForm &a_ho, const Array<int> &ess_tdof_list,
206  int ref_type=BasisType::GaussLobatto)
207  {
208  lor = new LORDiscretization(a_ho, ess_tdof_list, ref_type);
210  }
211 
212 #ifdef MFEM_USE_MPI
213  /// @brief Create a solver of type @a SolverType, formed using the assembled
214  /// HypreParMatrix of the LOR version of @a a_ho. @see ParLORDiscretization
215  LORSolver(ParBilinearForm &a_ho, const Array<int> &ess_tdof_list,
216  int ref_type=BasisType::GaussLobatto)
217  {
218  lor = new ParLORDiscretization(a_ho, ess_tdof_list, ref_type);
220  }
221 #endif
222 
223  /// @brief Create a solver of type @a SolverType using Operator @a op and
224  /// arguments @a args.
225  template <typename... Args>
226  LORSolver(const Operator &op, LORBase &lor_, Args&&... args) : solver(args...)
227  {
228  lor = &lor_;
229  own_lor = false;
230  SetOperator(op);
231  }
232 
233  /// @brief Create a solver of type @a SolverType using the assembled LOR
234  /// operator represented by @a lor_.
235  ///
236  /// The given @a args will be used as arguments to the solver constructor.
237  template <typename... Args>
238  LORSolver(LORBase &lor_, Args&&... args)
239  : LORSolver(*lor_.GetAssembledSystem(), lor_, args...) { }
240 
241  void SetOperator(const Operator &op)
242  {
243  solver.SetOperator(op);
244  width = solver.Width();
245  height = solver.Height();
246  }
247 
248  void Mult(const Vector &x, Vector &y) const { solver.Mult(x, y); }
249 
250  /// Access the underlying solver.
251  SolverType &GetSolver() { return solver; }
252 
253  /// Access the underlying solver.
254  const SolverType &GetSolver() const { return solver; }
255 
256  /// Access the LOR discretization object.
257  const LORBase &GetLOR() const { return *lor; }
258 
259  ~LORSolver() { if (own_lor) { delete lor; } }
260 };
261 
262 } // namespace mfem
263 
264 #endif
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Definition: lor.hpp:22
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...
Definition: lor.cpp:376
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system.
Definition: lor.cpp:258
Create and assemble a low-order refined version of a ParBilinearForm.
Definition: lor.hpp:156
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
FESpaceType GetFESpaceType() const
Returns the type of finite element space: H1, ND, RT or L2.
Definition: lor.cpp:68
LORBase * lor
Definition: lor.hpp:198
int GetLOROrder() const
Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
Definition: lor.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
FiniteElementSpace & fes_ho
Definition: lor.hpp:65
Container class for integration rules.
Definition: intrules.hpp:311
const Array< int > & GetDofPermutation() const
Returns the permutation that maps LOR DOFs to high-order DOFs.
Definition: lor.cpp:246
Vector px
Definition: lor.hpp:201
Abstract parallel finite element space.
Definition: pfespace.hpp:28
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Definition: lor.cpp:462
Create and assemble a low-order refined version of a BilinearForm.
Definition: lor.hpp:127
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...
Definition: lor.hpp:215
bool HasSameDofNumbering() const
Definition: lor.cpp:252
void AssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho.
Definition: lor.cpp:407
FiniteElementSpace * fes
Definition: lor.hpp:68
Data type sparse matrix.
Definition: sparsemat.hpp:46
Mesh * mesh
Definition: lor.hpp:66
FiniteElementSpace & GetFESpace() const
Returns the low-order refined finite element space.
Definition: lor.hpp:121
void SetupProlongationAndRestriction()
Definition: lor.cpp:287
SolverType solver
Definition: lor.hpp:200
Array< int > perm
Definition: lor.hpp:71
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: lor.hpp:241
LORSolver(LORBase &lor_, Args &&...args)
Create a solver of type SolverType using the assembled LOR operator represented by lor_...
Definition: lor.hpp:238
BilinearForm * a
Definition: lor.hpp:69
SolverType & GetSolver()
Access the underlying solver.
Definition: lor.hpp:251
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
bool own_lor
Definition: lor.hpp:199
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
LORSolver(const Operator &op, LORBase &lor_, Args &&...args)
Create a solver of type SolverType using Operator op and arguments args.
Definition: lor.hpp:226
void AssembleSystem(ParBilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho.
Definition: lor.cpp:454
const SolverType & GetSolver() const
Access the underlying solver.
Definition: lor.hpp:254
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: lor.hpp:248
void ConstructDofPermutation() const
Definition: lor.cpp:204
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
ParFiniteElementSpace & GetParFESpace() const
Return the LOR ParFiniteElementSpace.
Definition: lor.cpp:468
FiniteElementCollection * fec
Definition: lor.hpp:67
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...
Definition: lor.cpp:423
void ConstructLocalDofPermutation(Array< int > &perm_) const
Definition: lor.cpp:85
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...
Definition: lor.hpp:205
Class for parallel bilinear form.
LORBase(FiniteElementSpace &fes_ho_)
Definition: lor.cpp:348
SparseMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a SparseMatrix.
Definition: lor.cpp:415
Vector data type.
Definition: vector.hpp:60
Vector py
Definition: lor.hpp:201
void AssembleSystem_(BilinearForm &a_ho, const Array< int > &ess_dofs)
Definition: lor.cpp:264
Base class for solvers.
Definition: operator.hpp:651
const LORBase & GetLOR() const
Access the LOR discretization object.
Definition: lor.hpp:257
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
OperatorHandle A
Definition: lor.hpp:70
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
Represents a solver of type SolverType created using the low-order refined version of the given Bilin...
Definition: lor.hpp:195