MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilinearform_ext.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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_BILINEARFORM_EXT
13#define MFEM_BILINEARFORM_EXT
14
15#include "../config/config.hpp"
16#include "fespace.hpp"
17#include "../general/device.hpp"
18
19namespace mfem
20{
21
22class BilinearForm;
23class MixedBilinearForm;
24class DiscreteLinearOperator;
25
26/// Class extending the BilinearForm class to support different AssemblyLevels.
27/** FA - Full Assembly
28 PA - Partial Assembly
29 EA - Element Assembly
30 MF - Matrix Free
31*/
33{
34protected:
35 BilinearForm *a; ///< Not owned
36
37public:
39
40 MemoryClass GetMemoryClass() const override
42
43 /// Get the finite element space prolongation matrix
44 const Operator *GetProlongation() const override;
45
46 /// Get the finite element space restriction matrix
47 const Operator *GetRestriction() const override;
48
49 /// Assemble at the level given for the BilinearFormExtension subclass
50 virtual void Assemble() = 0;
51
52 void AssembleDiagonal(Vector &diag) const override
53 {
54 MFEM_ABORT("AssembleDiagonal not implemented for this assembly level!");
55 }
56
57 virtual void FormSystemMatrix(const Array<int> &ess_tdof_list,
58 OperatorHandle &A) = 0;
59 virtual void FormLinearSystem(const Array<int> &ess_tdof_list,
60 Vector &x, Vector &b,
61 OperatorHandle &A, Vector &X, Vector &B,
62 int copy_interior = 0) = 0;
63 virtual void Update() = 0;
64};
65
66/// Data and methods for partially-assembled bilinear forms
68{
69protected:
70 const FiniteElementSpace *trial_fes, *test_fes; // Not owned
71 /// Attributes of all mesh elements.
73 mutable Vector tmp_evec; // Work array
79 const Operator *elem_restrict; // Not owned
82
83public:
85
86 void Assemble() override;
87 void AssembleDiagonal(Vector &diag) const override;
88 void FormSystemMatrix(const Array<int> &ess_tdof_list,
89 OperatorHandle &A) override;
90 void FormLinearSystem(const Array<int> &ess_tdof_list,
91 Vector &x, Vector &b,
92 OperatorHandle &A, Vector &X, Vector &B,
93 int copy_interior = 0) override;
94 void Mult(const Vector &x, Vector &y) const override;
95 void MultTranspose(const Vector &x, Vector &y) const override;
96 void Update() override;
97
98protected:
100
101 /// @brief Accumulate the action (or transpose) of the integrator on @a x
102 /// into @a y, taking into account the (possibly null) @a markers array.
103 ///
104 /// If @a markers is non-null, then only those elements or boundary elements
105 /// whose attribute is marked in the markers array will be added to @a y.
106 ///
107 /// @param integ The integrator (domain, boundary, or boundary face).
108 /// @param x Input E-vector.
109 /// @param markers Marked attributes (possibly null, meaning all attributes).
110 /// @param attributes Array of element or boundary element attributes.
111 /// @param transpose Compute the action or transpose of the integrator .
112 /// @param y Output E-vector
114 const Vector &x,
115 const Array<int> *markers,
116 const Array<int> &attributes,
117 const bool transpose,
118 Vector &y) const;
119
120 /// @brief Performs the same function as AddMultWithMarkers, but takes as
121 /// input and output face normal derivatives.
122 ///
123 /// This is required when the integrator requires face normal derivatives,
124 /// for example, DGDiffusionIntegrator.
125 ///
126 /// This is called when the integrator's member function
127 /// BilinearFormIntegrator::RequiresFaceNormalDerivatives() returns true.
129 const BilinearFormIntegrator &integ,
130 const Vector &x,
131 const Vector &dxdn,
132 const Array<int> *markers,
133 const Array<int> &attributes,
134 Vector &y,
135 Vector &dydn) const;
136};
137
138/// Data and methods for element-assembled bilinear forms
140{
141protected:
142 int ne;
144 // The element matrices are stored row major
150
151public:
153
154 void Assemble() override;
155 void Mult(const Vector &x, Vector &y) const override;
156 void MultTranspose(const Vector &x, Vector &y) const override;
157
158 /// @brief Populates @a element_matrices with the element matrices.
159 ///
160 /// The element matrices are converted from row-major (how they are stored in
161 /// @a ea_data) to column-major format.
162 ///
163 /// If @a ordering is ElementDofOrdering::NATIVE, then the matrices are
164 /// reordered from the lexicographic ordering used internally.
165 void GetElementMatrices(DenseTensor &element_matrices,
166 ElementDofOrdering ordering,
167 bool add_bdr);
168};
169
170/// Data and methods for fully-assembled bilinear forms
172{
173private:
174 SparseMatrix *mat;
175 mutable Vector dg_x, dg_y;
176
177public:
179
180 void Assemble() override;
181 void RAP(OperatorHandle &A);
182 /** @note Always does `DIAG_ONE` policy to be consistent with
183 `Operator::FormConstrainedSystemOperator`. */
184 void EliminateBC(const Array<int> &ess_dofs, OperatorHandle &A);
185 void FormSystemMatrix(const Array<int> &ess_tdof_list,
186 OperatorHandle &A) override;
187 void FormLinearSystem(const Array<int> &ess_tdof_list,
188 Vector &x, Vector &b,
189 OperatorHandle &A, Vector &X, Vector &B,
190 int copy_interior = 0) override;
191 void Mult(const Vector &x, Vector &y) const override;
192 void MultTranspose(const Vector &x, Vector &y) const override;
193
194 /** DGMult and DGMultTranspose use the extended L-vector to perform the
195 computation. */
196 void DGMult(const Vector &x, Vector &y) const;
197 void DGMultTranspose(const Vector &x, Vector &y) const;
198};
199
200/// Data and methods for matrix-free bilinear forms
202{
203protected:
204 const FiniteElementSpace *trial_fes, *test_fes; // Not owned
208 const Operator *elem_restrict; // Not owned
211
212public:
214
215 void Assemble() override;
216 void AssembleDiagonal(Vector &diag) const override;
217 void FormSystemMatrix(const Array<int> &ess_tdof_list,
218 OperatorHandle &A) override;
219 void FormLinearSystem(const Array<int> &ess_tdof_list,
220 Vector &x, Vector &b,
221 OperatorHandle &A, Vector &X, Vector &B,
222 int copy_interior = 0) override;
223 void Mult(const Vector &x, Vector &y) const override;
224 void MultTranspose(const Vector &x, Vector &y) const override;
225 void Update() override;
226};
227
228/// Class extending the MixedBilinearForm class to support different AssemblyLevels.
229/** FA - Full Assembly
230 PA - Partial Assembly
231 EA - Element Assembly
232 MF - Matrix Free
233*/
235{
236protected:
237 MixedBilinearForm *a; ///< Not owned
238
239public:
241
243 { return Device::GetMemoryClass(); }
244
245 /// Get the finite element space prolongation matrix
246 const Operator *GetProlongation() const override;
247
248 /// Get the finite element space restriction matrix
249 const Operator *GetRestriction() const override;
250
251 /// Get the output finite element space restriction matrix
252 const Operator *GetOutputProlongation() const override;
253
254 /// Get the output finite element space restriction matrix
255 const Operator *GetOutputRestriction() const override;
256
257 virtual void Assemble() = 0;
258 virtual void FormRectangularSystemOperator(const Array<int> &trial_tdof_list,
259 const Array<int> &test_tdof_list,
260 OperatorHandle &A) = 0;
261 virtual void FormRectangularLinearSystem(const Array<int> &trial_tdof_list,
262 const Array<int> &test_tdof_list,
263 Vector &x, Vector &b,
264 OperatorHandle &A, Vector &X, Vector &B) = 0;
265
266 virtual void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const = 0;
267
268 virtual void Update() = 0;
269};
270
271/// Data and methods for partially-assembled mixed bilinear forms
273{
274protected:
275 const FiniteElementSpace *trial_fes, *test_fes; // Not owned
277 const Operator *elem_restrict_trial; // Not owned
278 const Operator *elem_restrict_test; // Not owned
279
280 /// Helper function to set up inputs/outputs for Mult or MultTranspose
281 void SetupMultInputs(const Operator *elem_restrict_x,
282 const Vector &x, Vector &localX,
283 const Operator *elem_restrict_y,
284 Vector &y, Vector &localY, const real_t c) const;
285
286public:
288
289 /// Partial assembly of all internal integrators
290 void Assemble() override;
291 /**
292 @brief Setup OperatorHandle A to contain constrained linear operator
293
294 OperatorHandle A contains matrix-free constrained operator formed for RAP
295 system where ess_tdof_list are in trial space and eliminated from
296 "columns" of A.
297 */
298 void FormRectangularSystemOperator(const Array<int> &trial_tdof_list,
299 const Array<int> &test_tdof_list,
300 OperatorHandle &A) override;
301 /**
302 Setup OperatorHandle A to contain constrained linear operator and
303 eliminate columns corresponding to essential dofs from system,
304 updating RHS B vector with the results.
305 */
306 void FormRectangularLinearSystem(const Array<int> &trial_tdof_list,
307 const Array<int> &test_tdof_list,
308 Vector &x, Vector &b,
309 OperatorHandle &A, Vector &X, Vector &B) override;
310 /// y = A*x
311 void Mult(const Vector &x, Vector &y) const override;
312 /// y += c*A*x
313 void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const override;
314 /// y = A^T*x
315 void MultTranspose(const Vector &x, Vector &y) const override;
316 /// y += c*A^T*x
317 void AddMultTranspose(const Vector &x, Vector &y,
318 const real_t c=1.0) const override;
319 /// Assemble the diagonal of ADA^T for a diagonal vector D.
320 void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const override;
321
322 /// Update internals for when a new MixedBilinearForm is given to this class
323 void Update() override;
324};
325
326
327/**
328 @brief Partial assembly extension for DiscreteLinearOperator
329
330 This acts very much like PAMixedBilinearFormExtension, but its
331 FormRectangularSystemOperator implementation emulates 'Set' rather than
332 'Add' in the assembly case.
333*/
335{
336public:
338
339 /// Partial assembly of all internal integrators
340 void Assemble() override;
341
342 void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const override;
343
344 void AddMultTranspose(const Vector &x, Vector &y,
345 const real_t c=1.0) const override;
346
348 OperatorHandle& A) override;
349
350 const Operator * GetOutputRestrictionTranspose() const override;
351
352private:
353 Vector test_multiplicity;
354};
355
356}
357
358#endif
Class extending the BilinearForm class to support different AssemblyLevels.
virtual void Assemble()=0
Assemble at the level given for the BilinearFormExtension subclass.
BilinearFormExtension(BilinearForm *form)
void AssembleDiagonal(Vector &diag) const override
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)=0
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)=0
const Operator * GetRestriction() const override
Get the finite element space restriction matrix.
MemoryClass GetMemoryClass() const override
Return the MemoryClass preferred by the Operator.
BilinearForm * a
Not owned.
const Operator * GetProlongation() const override
Get the finite element space prolongation matrix.
Abstract base class BilinearFormIntegrator.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Rank 3 tensor (array of matrices)
static MemoryClass GetMemoryClass()
(DEPRECATED) Equivalent to GetDeviceMemoryClass().
Definition device.hpp:286
static MemoryClass GetDeviceMemoryClass()
Get the current Device MemoryClass. This is the MemoryClass used by most MFEM device kernels to acces...
Definition device.hpp:282
Data and methods for element-assembled bilinear forms.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
void GetElementMatrices(DenseTensor &element_matrices, ElementDofOrdering ordering, bool add_bdr)
Populates element_matrices with the element matrices.
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
EABilinearFormExtension(BilinearForm *form)
void Assemble() override
Assemble at the level given for the BilinearFormExtension subclass.
Data and methods for fully-assembled bilinear forms.
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void DGMultTranspose(const Vector &x, Vector &y) const
void EliminateBC(const Array< int > &ess_dofs, OperatorHandle &A)
void DGMult(const Vector &x, Vector &y) const
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A) override
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0) override
void Assemble() override
Assemble at the level given for the BilinearFormExtension subclass.
FABilinearFormExtension(BilinearForm *form)
Base class for operators that extracts Face degrees of freedom.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
Data and methods for matrix-free bilinear forms.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
const FiniteElementSpace * trial_fes
void AssembleDiagonal(Vector &diag) const override
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
const FaceRestriction * bdr_face_restrict_lex
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0) override
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void Assemble() override
Assemble at the level given for the BilinearFormExtension subclass.
const FiniteElementSpace * test_fes
const FaceRestriction * int_face_restrict_lex
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A) override
MFBilinearFormExtension(BilinearForm *form)
Class extending the MixedBilinearForm class to support different AssemblyLevels.
MixedBilinearForm * a
Not owned.
virtual void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const =0
const Operator * GetOutputProlongation() const override
Get the output finite element space restriction matrix.
virtual void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)=0
const Operator * GetProlongation() const override
Get the finite element space prolongation matrix.
MemoryClass GetMemoryClass() const override
Return the MemoryClass preferred by the Operator.
const Operator * GetOutputRestriction() const override
Get the output finite element space restriction matrix.
MixedBilinearFormExtension(MixedBilinearForm *form)
virtual void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)=0
const Operator * GetRestriction() const override
Get the finite element space restriction matrix.
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Abstract operator.
Definition operator.hpp:25
Data and methods for partially-assembled bilinear forms.
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A) override
Array< int > elem_attributes
Attributes of all mesh elements.
void AssembleDiagonal(Vector &diag) const override
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
const FiniteElementSpace * trial_fes
const FaceRestriction * int_face_restrict_lex
const FaceRestriction * bdr_face_restrict_lex
void AddMultNormalDerivativesWithMarkers(const BilinearFormIntegrator &integ, const Vector &x, const Vector &dxdn, const Array< int > *markers, const Array< int > &attributes, Vector &y, Vector &dydn) const
Performs the same function as AddMultWithMarkers, but takes as input and output face normal derivativ...
void AddMultWithMarkers(const BilinearFormIntegrator &integ, const Vector &x, const Array< int > *markers, const Array< int > &attributes, const bool transpose, Vector &y) const
Accumulate the action (or transpose) of the integrator on x into y, taking into account the (possibly...
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0) override
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
const FiniteElementSpace * test_fes
void Assemble() override
Assemble at the level given for the BilinearFormExtension subclass.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
void SetupRestrictionOperators(const L2FaceValues m)
Partial assembly extension for DiscreteLinearOperator.
void Assemble() override
Partial assembly of all internal integrators.
PADiscreteLinearOperatorExtension(DiscreteLinearOperator *linop)
void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const override
Operator application: y+=A(x) (default) or y+=a*A(x).
void AddMultTranspose(const Vector &x, Vector &y, const real_t c=1.0) const override
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
const Operator * GetOutputRestrictionTranspose() const override
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type...
void FormRectangularSystemOperator(const Array< int > &, const Array< int > &, OperatorHandle &A) override
Data and methods for partially-assembled mixed bilinear forms.
void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A) override
Setup OperatorHandle A to contain constrained linear operator.
const FiniteElementSpace * test_fes
void Assemble() override
Partial assembly of all internal integrators.
PAMixedBilinearFormExtension(MixedBilinearForm *form)
void MultTranspose(const Vector &x, Vector &y) const override
y = A^T*x
const FiniteElementSpace * trial_fes
void AddMultTranspose(const Vector &x, Vector &y, const real_t c=1.0) const override
y += c*A^T*x
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const override
Assemble the diagonal of ADA^T for a diagonal vector D.
void SetupMultInputs(const Operator *elem_restrict_x, const Vector &x, Vector &localX, const Operator *elem_restrict_y, Vector &y, Vector &localY, const real_t c) const
Helper function to set up inputs/outputs for Mult or MultTranspose.
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B) override
void Update() override
Update internals for when a new MixedBilinearForm is given to this class.
void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const override
y += c*A*x
void Mult(const Vector &x, Vector &y) const override
y = A*x
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:82
real_t b
Definition lissajous.cpp:42
MemoryClass
Memory classes identify sets of memory types.
float real_t
Definition config.hpp:43
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:83