MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
bilinearform_ext.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_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
42
43 /// Get the finite element space prolongation matrix
44 virtual const Operator *GetProlongation() const;
45
46 /// Get the finite element space restriction matrix
47 virtual const Operator *GetRestriction() const;
48
49 /// Assemble at the level given for the BilinearFormExtension subclass
50 virtual void Assemble() = 0;
51
52 virtual void AssembleDiagonal(Vector &diag) const
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();
87 void AssembleDiagonal(Vector &diag) const;
88 void FormSystemMatrix(const Array<int> &ess_tdof_list, OperatorHandle &A);
89 void FormLinearSystem(const Array<int> &ess_tdof_list,
90 Vector &x, Vector &b,
91 OperatorHandle &A, Vector &X, Vector &B,
92 int copy_interior = 0);
93 void Mult(const Vector &x, Vector &y) const;
94 void MultTranspose(const Vector &x, Vector &y) const;
95 void Update();
96
97protected:
99
100 /// @brief Accumulate the action (or transpose) of the integrator on @a x
101 /// into @a y, taking into account the (possibly null) @a markers array.
102 ///
103 /// If @a markers is non-null, then only those elements or boundary elements
104 /// whose attribute is marked in the markers array will be added to @a y.
105 ///
106 /// @param integ The integrator (domain, boundary, or boundary face).
107 /// @param x Input E-vector.
108 /// @param markers Marked attributes (possibly null, meaning all attributes).
109 /// @param attributes Array of element or boundary element attributes.
110 /// @param transpose Compute the action or transpose of the integrator .
111 /// @param y Output E-vector
113 const Vector &x,
114 const Array<int> *markers,
115 const Array<int> &attributes,
116 const bool transpose,
117 Vector &y) const;
118
119 /// @brief Performs the same function as AddMultWithMarkers, but takes as
120 /// input and output face normal derivatives.
121 ///
122 /// This is required when the integrator requires face normal derivatives,
123 /// for example, DGDiffusionIntegrator.
124 ///
125 /// This is called when the integrator's member function
126 /// BilinearFormIntegrator::RequiresFaceNormalDerivatives() returns true.
128 const BilinearFormIntegrator &integ,
129 const Vector &x,
130 const Vector &dxdn,
131 const Array<int> *markers,
132 const Array<int> &attributes,
133 Vector &y,
134 Vector &dydn) const;
135};
136
137/// Data and methods for element-assembled bilinear forms
139{
140protected:
141 int ne;
143 // The element matrices are stored row major
149
150public:
152
153 void Assemble();
154 void Mult(const Vector &x, Vector &y) const;
155 void MultTranspose(const Vector &x, Vector &y) const;
156};
157
158/// Data and methods for fully-assembled bilinear forms
160{
161private:
162 SparseMatrix *mat;
163 mutable Vector dg_x, dg_y;
164
165public:
167
168 void Assemble();
169 void RAP(OperatorHandle &A);
170 /** @note Always does `DIAG_ONE` policy to be consistent with
171 `Operator::FormConstrainedSystemOperator`. */
172 void EliminateBC(const Array<int> &ess_dofs, OperatorHandle &A);
173 void FormSystemMatrix(const Array<int> &ess_tdof_list, OperatorHandle &A);
174 void FormLinearSystem(const Array<int> &ess_tdof_list,
175 Vector &x, Vector &b,
176 OperatorHandle &A, Vector &X, Vector &B,
177 int copy_interior = 0);
178 void Mult(const Vector &x, Vector &y) const;
179 void MultTranspose(const Vector &x, Vector &y) const;
180
181 /** DGMult and DGMultTranspose use the extended L-vector to perform the
182 computation. */
183 void DGMult(const Vector &x, Vector &y) const;
184 void DGMultTranspose(const Vector &x, Vector &y) const;
185};
186
187/// Data and methods for matrix-free bilinear forms
189{
190protected:
191 const FiniteElementSpace *trial_fes, *test_fes; // Not owned
195 const Operator *elem_restrict; // Not owned
198
199public:
201
202 void Assemble();
203 void AssembleDiagonal(Vector &diag) const;
204 void FormSystemMatrix(const Array<int> &ess_tdof_list, OperatorHandle &A);
205 void FormLinearSystem(const Array<int> &ess_tdof_list,
206 Vector &x, Vector &b,
207 OperatorHandle &A, Vector &X, Vector &B,
208 int copy_interior = 0);
209 void Mult(const Vector &x, Vector &y) const;
210 void MultTranspose(const Vector &x, Vector &y) const;
211 void Update();
212};
213
214/// Class extending the MixedBilinearForm class to support different AssemblyLevels.
215/** FA - Full Assembly
216 PA - Partial Assembly
217 EA - Element Assembly
218 MF - Matrix Free
219*/
221{
222protected:
223 MixedBilinearForm *a; ///< Not owned
224
225public:
227
229 { return Device::GetMemoryClass(); }
230
231 /// Get the finite element space prolongation matrix
232 virtual const Operator *GetProlongation() const;
233
234 /// Get the finite element space restriction matrix
235 virtual const Operator *GetRestriction() const;
236
237 /// Get the output finite element space restriction matrix
238 virtual const Operator *GetOutputProlongation() const;
239
240 /// Get the output finite element space restriction matrix
241 virtual const Operator *GetOutputRestriction() const;
242
243 virtual void Assemble() = 0;
244 virtual void FormRectangularSystemOperator(const Array<int> &trial_tdof_list,
245 const Array<int> &test_tdof_list,
246 OperatorHandle &A) = 0;
247 virtual void FormRectangularLinearSystem(const Array<int> &trial_tdof_list,
248 const Array<int> &test_tdof_list,
249 Vector &x, Vector &b,
250 OperatorHandle &A, Vector &X, Vector &B) = 0;
251
252 virtual void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const = 0;
253
254 virtual void Update() = 0;
255};
256
257/// Data and methods for partially-assembled mixed bilinear forms
259{
260protected:
261 const FiniteElementSpace *trial_fes, *test_fes; // Not owned
263 const Operator *elem_restrict_trial; // Not owned
264 const Operator *elem_restrict_test; // Not owned
265
266 /// Helper function to set up inputs/outputs for Mult or MultTranspose
267 void SetupMultInputs(const Operator *elem_restrict_x,
268 const Vector &x, Vector &localX,
269 const Operator *elem_restrict_y,
270 Vector &y, Vector &localY, const real_t c) const;
271
272public:
274
275 /// Partial assembly of all internal integrators
276 void Assemble();
277 /**
278 @brief Setup OperatorHandle A to contain constrained linear operator
279
280 OperatorHandle A contains matrix-free constrained operator formed for RAP
281 system where ess_tdof_list are in trial space and eliminated from
282 "columns" of A.
283 */
284 void FormRectangularSystemOperator(const Array<int> &trial_tdof_list,
285 const Array<int> &test_tdof_list,
286 OperatorHandle &A);
287 /**
288 Setup OperatorHandle A to contain constrained linear operator and
289 eliminate columns corresponding to essential dofs from system,
290 updating RHS B vector with the results.
291 */
292 void FormRectangularLinearSystem(const Array<int> &trial_tdof_list,
293 const Array<int> &test_tdof_list,
294 Vector &x, Vector &b,
295 OperatorHandle &A, Vector &X, Vector &B);
296 /// y = A*x
297 void Mult(const Vector &x, Vector &y) const;
298 /// y += c*A*x
299 void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const;
300 /// y = A^T*x
301 void MultTranspose(const Vector &x, Vector &y) const;
302 /// y += c*A^T*x
303 void AddMultTranspose(const Vector &x, Vector &y, const real_t c=1.0) const;
304 /// Assemble the diagonal of ADA^T for a diagonal vector D.
305 void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const;
306
307 /// Update internals for when a new MixedBilinearForm is given to this class
308 void Update();
309};
310
311
312/**
313 @brief Partial assembly extension for DiscreteLinearOperator
314
315 This acts very much like PAMixedBilinearFormExtension, but its
316 FormRectangularSystemOperator implementation emulates 'Set' rather than
317 'Add' in the assembly case.
318*/
320{
321public:
323
324 /// Partial assembly of all internal integrators
325 void Assemble();
326
327 void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const;
328
329 void AddMultTranspose(const Vector &x, Vector &y, const real_t c=1.0) const;
330
332 OperatorHandle& A);
333
335
336private:
337 Vector test_multiplicity;
338};
339
340}
341
342#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)
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
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
BilinearForm * a
Not owned.
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Abstract base class BilinearFormIntegrator.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
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 Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
EABilinearFormExtension(BilinearForm *form)
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Data and methods for fully-assembled bilinear forms.
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
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 MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
FABilinearFormExtension(BilinearForm *form)
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
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:220
Data and methods for matrix-free bilinear forms.
void Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
const FiniteElementSpace * trial_fes
const FaceRestriction * bdr_face_restrict_lex
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
const FiniteElementSpace * test_fes
const FaceRestriction * int_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)
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
MFBilinearFormExtension(BilinearForm *form)
void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
Class extending the MixedBilinearForm class to support different AssemblyLevels.
MixedBilinearForm * a
Not owned.
virtual const Operator * GetOutputProlongation() const
Get the output finite element space restriction matrix.
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
virtual void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const =0
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
virtual void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)=0
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
virtual const Operator * GetOutputRestriction() const
Get the output 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 Assemble()
Assemble at the level given for the BilinearFormExtension subclass.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Array< int > elem_attributes
Attributes of all mesh elements.
const FiniteElementSpace * trial_fes
const FaceRestriction * int_face_restrict_lex
const FaceRestriction * bdr_face_restrict_lex
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
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 AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
const FiniteElementSpace * test_fes
void SetupRestrictionOperators(const L2FaceValues m)
Partial assembly extension for DiscreteLinearOperator.
PADiscreteLinearOperatorExtension(DiscreteLinearOperator *linop)
const Operator * GetOutputRestrictionTranspose() const
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type...
void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const
Operator application: y+=A(x) (default) or y+=a*A(x).
void FormRectangularSystemOperator(const Array< int > &, const Array< int > &, OperatorHandle &A)
void Assemble()
Partial assembly of all internal integrators.
void AddMultTranspose(const Vector &x, Vector &y, const real_t c=1.0) const
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
Data and methods for partially-assembled mixed bilinear forms.
const FiniteElementSpace * test_fes
void Assemble()
Partial assembly of all internal integrators.
PAMixedBilinearFormExtension(MixedBilinearForm *form)
void AddMultTranspose(const Vector &x, Vector &y, const real_t c=1.0) const
y += c*A^T*x
const FiniteElementSpace * trial_fes
void Mult(const Vector &x, Vector &y) const
y = A*x
void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Setup OperatorHandle A to contain constrained linear operator.
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 AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T for a diagonal vector D.
void AddMult(const Vector &x, Vector &y, const real_t c=1.0) const
y += c*A*x
void MultTranspose(const Vector &x, Vector &y) const
y = A^T*x
void Update()
Update internals for when a new MixedBilinearForm is given to this class.
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
real_t b
Definition lissajous.cpp:42
MemoryClass
Memory classes identify sets of memory types.
float real_t
Definition config.hpp:43