MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
lor.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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
17namespace mfem
18{
19
20/// @brief Abstract base class for LORDiscretization and ParLORDiscretization
21/// classes, which construct low-order refined versions of bilinear forms.
23{
24private:
25 using GetIntegratorsFn = Array<BilinearFormIntegrator*> *(BilinearForm::*)();
26 using GetMarkersFn = Array<Array<int>*> *(BilinearForm::*)();
27 using AddIntegratorFn = void (BilinearForm::*)(BilinearFormIntegrator*);
28 using AddIntegratorMarkersFn =
30
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
62protected:
63 enum FESpaceType { H1, ND, RT, L2, INVALID };
64
67 Mesh *mesh = nullptr;
70 BilinearForm *a = nullptr;
74
75 /// Constructs the local DOF (ldof) permutation. In parallel this is used as
76 /// an intermediate step in computing the DOF permutation (see
77 /// ConstructDofPermutation and GetDofPermutation).
79
80 /// Construct the permutation that maps LOR DOFs to high-order DOFs. See
81 /// GetDofPermutation.
82 void ConstructDofPermutation() const;
83
84 /// Returns true if the LOR space and HO space have the same DOF numbering
85 /// (H1 or L2 spaces), false otherwise (ND or RT spaces).
86 bool HasSameDofNumbering() const;
87
88 /// Sets up the prolongation and restriction operators required in the case
89 /// of different DOF numberings (ND or RT spaces) or nonconforming spaces.
91
92 /// Returns the type of finite element space: H1, ND, RT or L2.
94
95 /// Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
96 int GetLOROrder() const;
97
98 /// Construct the LOR space (overridden for serial and parallel versions).
99 virtual void FormLORSpace() = 0;
100
101 /// Construct the LORBase object for the given FE space and refinement type.
102 LORBase(FiniteElementSpace &fes_ho_, int ref_type_);
103
104public:
105 /// Returns the assembled LOR system (const version).
106 const OperatorHandle &GetAssembledSystem() const;
107
108 /// Returns the assembled LOR system (non-const version).
110
111 /// Assembles the LOR system corresponding to @a a_ho.
112 void AssembleSystem(BilinearForm &a_ho, const Array<int> &ess_dofs);
113
114 /// Assembles the LOR system corresponding to @a a_ho using the legacy method.
115 void LegacyAssembleSystem(BilinearForm &a_ho, const Array<int> &ess_dofs);
116
117 /// @brief Returns the permutation that maps LOR DOFs to high-order DOFs.
118 ///
119 /// This permutation is constructed the first time it is requested, and then
120 /// is cached. For H1 and L2 finite element spaces (or for nonconforming
121 /// spaces) this is the identity. In these cases, RequiresDofPermutation will
122 /// return false. However, if the DOF permutation is requested, an identity
123 /// permutation will be built and returned.
124 ///
125 /// For vector finite element spaces (ND and RT), the DOF permutation is
126 /// nontrivial. Returns an array @a perm such that, given an index @a i of a
127 /// LOR dof, @a perm[i] is the index of the corresponding HO dof.
128 const Array<int> &GetDofPermutation() const;
129
130 /// Returns the low-order refined finite element space.
132
133 virtual ~LORBase();
134};
135
136/// Create and assemble a low-order refined version of a BilinearForm.
138{
139protected:
140 void FormLORSpace() override;
141public:
142 /// @brief Construct the low-order refined version of @a a_ho using the given
143 /// list of essential DOFs.
144 ///
145 /// The mesh is refined using the refinement type specified by @a ref_type
146 /// (see Mesh::MakeRefined).
147 LORDiscretization(BilinearForm &a_ho, const Array<int> &ess_tdof_list,
149
150 /// @brief Construct a low-order refined version of the FiniteElementSpace @a
151 /// fes_ho.
152 ///
153 /// The mesh is refined using the refinement type specified by @a ref_type
154 /// (see Mesh::MakeRefined).
157
158 /// Return the assembled LOR operator as a SparseMatrix.
160};
161
162#ifdef MFEM_USE_MPI
163
164/// Create and assemble a low-order refined version of a ParBilinearForm.
166{
167protected:
168 void FormLORSpace() override;
169public:
170 /// @brief Construct the low-order refined version of @a a_ho using the given
171 /// list of essential DOFs.
172 ///
173 /// The mesh is refined using the refinement type specified by @a ref_type
174 /// (see ParMesh::MakeRefined).
175 ParLORDiscretization(ParBilinearForm &a_ho, const Array<int> &ess_tdof_list,
177
178 /// @brief Construct a low-order refined version of the ParFiniteElementSpace
179 /// @a pfes_ho.
180 ///
181 /// The mesh is refined using the refinement type specified by @a ref_type
182 /// (see ParMesh::MakeRefined).
185
186 /// Return the assembled LOR operator as a HypreParMatrix.
188
189 /// Return the LOR ParFiniteElementSpace.
191};
192
193#endif
194
195/// @brief Represents a solver of type @a SolverType created using the low-order
196/// refined version of the given BilinearForm or ParBilinearForm.
197///
198/// @note To achieve good solver performance, the high-order finite element
199/// space should use BasisType::GaussLobatto for H1 discretizations, and basis
200/// pair (BasisType::GaussLobatto, BasisType::IntegratedGLL) for Nedelec and
201/// Raviart-Thomas elements.
202template <typename SolverType>
203class LORSolver : public Solver
204{
205protected:
207 bool own_lor = true;
208 SolverType solver;
209public:
210 /// @brief Create a solver of type @a SolverType, formed using the assembled
211 /// SparseMatrix of the LOR version of @a a_ho. @see LORDiscretization
212 LORSolver(BilinearForm &a_ho, const Array<int> &ess_tdof_list,
213 int ref_type=BasisType::GaussLobatto)
214 {
215 lor = new LORDiscretization(a_ho, ess_tdof_list, ref_type);
217 }
218
219#ifdef MFEM_USE_MPI
220 /// @brief Create a solver of type @a SolverType, formed using the assembled
221 /// HypreParMatrix of the LOR version of @a a_ho. @see ParLORDiscretization
222 LORSolver(ParBilinearForm &a_ho, const Array<int> &ess_tdof_list,
223 int ref_type=BasisType::GaussLobatto)
224 {
225 lor = new ParLORDiscretization(a_ho, ess_tdof_list, ref_type);
227 }
228#endif
229
230 /// @brief Create a solver of type @a SolverType using Operator @a op and
231 /// arguments @a args.
232 template <typename... Args>
233 LORSolver(const Operator &op, LORBase &lor_, Args&&... args) : solver(args...)
234 {
235 lor = &lor_;
236 own_lor = false;
237 SetOperator(op);
238 }
239
240 /// @brief Create a solver of type @a SolverType using the assembled LOR
241 /// operator represented by @a lor_.
242 ///
243 /// The given @a args will be used as arguments to the solver constructor.
244 template <typename... Args>
245 LORSolver(LORBase &lor_, Args&&... args)
246 : LORSolver(*lor_.GetAssembledSystem(), lor_, args...) { }
247
248 void SetOperator(const Operator &op)
249 {
250 solver.SetOperator(op);
251 width = solver.Width();
252 height = solver.Height();
253 }
254
255 void Mult(const Vector &x, Vector &y) const { solver.Mult(x, y); }
256
257 /// Access the underlying solver.
258 SolverType &GetSolver() { return solver; }
259
260 /// Access the underlying solver.
261 const SolverType &GetSolver() const { return solver; }
262
263 /// Access the LOR discretization object.
264 const LORBase &GetLOR() const { return *lor; }
265
266 ~LORSolver() { if (own_lor) { delete lor; } }
267};
268
269#ifdef MFEM_USE_MPI
270
271// Template specialization for batched LOR AMS (implementation in lor_ams.cpp)
272template <>
273class LORSolver<HypreAMS> : public Solver
274{
275protected:
276 OperatorHandle A; ///< The assembled system matrix.
277 Vector *xyz = nullptr; ///< Data for vertex coordinate vectors.
278 HypreAMS *solver = nullptr; ///< The underlying AMS solver.
279public:
280 /// @brief Creates the AMS solvers for the given form and essential DOFs.
281 ///
282 /// Assembles the LOR matrices for the form @a a_ho and the associated
283 /// discrete gradient matrix and vertex coordinate vectors.
284 LORSolver(ParBilinearForm &a_ho, const Array<int> &ess_tdof_list,
285 int ref_type=BasisType::GaussLobatto);
286
287 /// Calls HypreAMS::SetOperator.
288 void SetOperator(const Operator &op);
289
290 /// Apply the action of the AMS preconditioner.
291 void Mult(const Vector &x, Vector &y) const;
292
293 /// Access the underlying solver.
295
296 /// Access the underlying solver (const version).
297 const HypreAMS &GetSolver() const;
298
299 ~LORSolver();
300};
301
302// Template specialization for batched LOR ADS (implementation in lor_ads.cpp)
303template <>
304class LORSolver<HypreADS> : public Solver
305{
306protected:
307 OperatorHandle A; ///< The assembled system matrix.
308 Vector *xyz = nullptr; ///< Data for vertex coordinate vectors.
309 HypreADS *solver = nullptr; ///< The underlying ADS solver.
310public:
311 /// @brief Creates the ADS solvers for the given form and essential DOFs.
312 ///
313 /// Assembles the LOR matrices for the form @a a_ho and the associated
314 /// discrete gradient matrix and vertex coordinate vectors.
315 LORSolver(ParBilinearForm &a_ho, const Array<int> &ess_tdof_list,
316 int ref_type=BasisType::GaussLobatto);
317
318 /// Calls HypreADS::SetOperator.
319 void SetOperator(const Operator &op);
320
321 /// Apply the action of the AMS preconditioner.
322 void Mult(const Vector &x, Vector &y) const;
323
324 /// Access the underlying solver.
326
327 /// Access the underlying solver (const version).
328 const HypreADS &GetSolver() const;
329
330 ~LORSolver();
331};
332
333#endif
334
335} // namespace mfem
336
337#endif
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
Efficient batched assembly of LOR discretizations on device.
Abstract base class BilinearFormIntegrator.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
The Auxiliary-space Divergence Solver in hypre.
Definition hypre.hpp:1922
The Auxiliary-space Maxwell Solver in hypre.
Definition hypre.hpp:1845
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Container class for integration rules.
Definition intrules.hpp:416
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Definition lor.hpp:23
LORBase(FiniteElementSpace &fes_ho_, int ref_type_)
Construct the LORBase object for the given FE space and refinement type.
Definition lor.cpp:337
int GetLOROrder() const
Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
Definition lor.cpp:84
void ConstructLocalDofPermutation(Array< int > &perm_) const
Definition lor.cpp:90
bool HasSameDofNumbering() const
Definition lor.cpp:258
void ConstructDofPermutation() const
Definition lor.cpp:209
FiniteElementSpace * fes
Definition lor.hpp:69
void AssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho.
Definition lor.cpp:366
FESpaceType GetFESpaceType() const
Returns the type of finite element space: H1, ND, RT or L2.
Definition lor.cpp:73
FiniteElementCollection * fec
Definition lor.hpp:68
void SetupProlongationAndRestriction()
Definition lor.cpp:276
class BatchedLORAssembly * batched_lor
Definition lor.hpp:71
FiniteElementSpace & fes_ho
Definition lor.hpp:66
OperatorHandle A
Definition lor.hpp:72
BilinearForm * a
Definition lor.hpp:70
virtual void FormLORSpace()=0
Construct the LOR space (overridden for serial and parallel versions).
Array< int > perm
Definition lor.hpp:73
int ref_type
Definition lor.hpp:65
const Array< int > & GetDofPermutation() const
Returns the permutation that maps LOR DOFs to high-order DOFs.
Definition lor.cpp:252
FiniteElementSpace & GetFESpace() const
Returns the low-order refined finite element space.
Definition lor.cpp:357
virtual ~LORBase()
Definition lor.cpp:430
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system (const version).
Definition lor.cpp:270
Mesh * mesh
Definition lor.hpp:67
void LegacyAssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho using the legacy method.
Definition lor.cpp:386
Create and assemble a low-order refined version of a BilinearForm.
Definition lor.hpp:138
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:439
void FormLORSpace() override
Construct the LOR space (overridden for serial and parallel versions).
Definition lor.cpp:456
SparseMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a SparseMatrix.
Definition lor.cpp:476
OperatorHandle A
The assembled system matrix.
Definition lor.hpp:307
OperatorHandle A
The assembled system matrix.
Definition lor.hpp:276
Represents a solver of type SolverType created using the low-order refined version of the given Bilin...
Definition lor.hpp:204
LORSolver(LORBase &lor_, Args &&... args)
Create a solver of type SolverType using the assembled LOR operator represented by lor_.
Definition lor.hpp:245
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition lor.hpp:255
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition lor.hpp:248
LORSolver(const Operator &op, LORBase &lor_, Args &&... args)
Create a solver of type SolverType using Operator op and arguments args.
Definition lor.hpp:233
const SolverType & GetSolver() const
Access the underlying solver.
Definition lor.hpp:261
SolverType & GetSolver()
Access the underlying solver.
Definition lor.hpp:258
LORBase * lor
Definition lor.hpp:206
const LORBase & GetLOR() const
Access the LOR discretization object.
Definition lor.hpp:264
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:222
SolverType solver
Definition lor.hpp:208
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:212
Mesh data type.
Definition mesh.hpp:56
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Abstract operator.
Definition operator.hpp:25
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
Create and assemble a low-order refined version of a ParBilinearForm.
Definition lor.hpp:166
ParFiniteElementSpace & GetParFESpace() const
Return the LOR ParFiniteElementSpace.
Definition lor.cpp:528
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:484
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Definition lor.cpp:522
void FormLORSpace() override
Construct the LOR space (overridden for serial and parallel versions).
Definition lor.cpp:501
Base class for solvers.
Definition operator.hpp:683
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80