MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
doftrans.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_DOFTRANSFORM
13#define MFEM_DOFTRANSFORM
14
15#include "../config/config.hpp"
16#include "../linalg/linalg.hpp"
17#include "intrules.hpp"
18
19namespace mfem
20{
21
22/** The StatelessDofTransformation class is an abstract base class for a family
23 of transformations that map local degrees of freedom (DoFs), contained
24 within individual elements, to global degrees of freedom, stored within
25 GridFunction objects.
26
27 In this context "stateless" means that the concrete classes derived from
28 StatelessDofTransformation do not store information about the relative
29 orientations of the faces with respect to their neighboring elements. In
30 other words there is no information specific to a particular element (aside
31 from the element type e.g. tetrahedron, wedge, or pyramid). The
32 StatelessDofTransformation provides access to the transformation operators
33 for specific relative face orientations. These are useful, for example, when
34 relating DoFs associated with distinct overlapping meshes such as parent and
35 sub-meshes.
36
37 These transformations are necessary to ensure that basis functions in
38 neighboring (or overlapping) elements align correctly. Closely related but
39 complementary transformations are required for the entries stored in
40 LinearForm and BilinearForm objects. The StatelessDofTransformation class
41 is designed to apply the action of both of these types of DoF
42 transformations.
43
44 Let the "primal transformation" be given by the operator T. This means that
45 given a local element vector v the data that must be placed into a
46 GridFunction object is v_t = T * v.
47
48 We also need the inverse of the primal transformation T^{-1} so that we can
49 recover the local element vector from data read out of a GridFunction
50 e.g. v = T^{-1} * v_t.
51
52 We need to preserve the action of our linear forms applied to primal
53 vectors. In other words, if f is the local vector computed by a linear
54 form then f * v = f_t * v_t (where "*" represents an inner product of
55 vectors). This requires that f_t = T^{-T} * f i.e. the "dual transform" is
56 given by the transpose of the inverse of the primal transformation.
57
58 For bilinear forms we require that v^T * A * v = v_t^T * A_t * v_t. This
59 implies that A_t = T^{-T} * A * T^{-1}. This can be accomplished by
60 performing dual transformations of the rows and columns of the matrix A.
61
62 For discrete linear operators the range must be modified with the primal
63 transformation rather than the dual transformation because the result is a
64 primal vector rather than a dual vector. This leads to the transformation
65 D_t = T * D * T^{-1}. This can be accomplished by using a primal
66 transformation on the columns of D and a dual transformation on its rows.
67*/
69{
70protected:
71 int size_;
72
74 : size_(size) {}
75
76public:
77 inline int Size() const { return size_; }
78 inline int Height() const { return size_; }
79 inline int NumRows() const { return size_; }
80 inline int Width() const { return size_; }
81 inline int NumCols() const { return size_; }
82
83 /// If the DofTransformation performs no transformation
84 virtual bool IsIdentity() const = 0;
85
86 /** Transform local DoFs to align with the global DoFs. For example, this
87 transformation can be used to map the local vector computed by
88 FiniteElement::Project() to the transformed vector stored within a
89 GridFunction object. */
90 virtual void TransformPrimal(const Array<int> & face_orientation,
91 real_t *v) const = 0;
92 inline void TransformPrimal(const Array<int> & face_orientation,
93 Vector &v) const
94 { TransformPrimal(face_orientation, v.GetData()); }
95
96 /** Inverse transform local DoFs. Used to transform DoFs from a global vector
97 back to their element-local form. For example, this must be used to
98 transform the vector obtained using GridFunction::GetSubVector before it
99 can be used to compute a local interpolation.
100 */
101 virtual void InvTransformPrimal(const Array<int> & face_orientation,
102 real_t *v) const = 0;
103 inline void InvTransformPrimal(const Array<int> & face_orientation,
104 Vector &v) const
105 { InvTransformPrimal(face_orientation, v.GetData()); }
106
107 /** Transform dual DoFs as computed by a LinearFormIntegrator before summing
108 into a LinearForm object. */
109 virtual void TransformDual(const Array<int> & face_orientation,
110 real_t *v) const = 0;
111 inline void TransformDual(const Array<int> & face_orientation,
112 Vector &v) const
113 { TransformDual(face_orientation, v.GetData()); }
114
115 /** Inverse Transform dual DoFs */
116 virtual void InvTransformDual(const Array<int> & face_orientation,
117 real_t *v) const = 0;
118 inline void InvTransformDual(const Array<int> & face_orientation,
119 Vector &v) const
120 { InvTransformDual(face_orientation, v.GetData()); }
121
122 virtual ~StatelessDofTransformation() = default;
123};
124
125/** The DofTransformation class is an extension of the
126 StatelessDofTransformation which stores the face orientations used to
127 select the necessary transformations which allows it to offer a collection
128 of convenience methods.
129
130 DofTransformation objects are provided by the FiniteElementSpace which has
131 access to the mesh and can therefore provide the face orientations. This is
132 convenient when working with GridFunction, LinearForm, or BilinearForm
133 objects or their parallel counterparts.
134
135 StatelessDofTransformation objects are provided by FiniteElement or
136 FiniteElementCollection objects which do not have access to face
137 orientation information. This can be useful in non-standard contexts such as
138 transferring finite element degrees of freedom between different meshes.
139 For examples of its use see the TransferMap used by the SubMesh class.
140 */
142{
143protected:
146 int vdim_;
148
149public:
150 /** @brief Default constructor which requires that SetDofTransformation be
151 called before use. */
152 DofTransformation(int vdim = 1, int ordering = 0)
153 : dof_trans_(NULL)
154 , vdim_(vdim)
155 , ordering_(ordering)
156 {}
157
158 /// Constructor with a known StatelessDofTransformation
160 int vdim = 1, int ordering = 0)
161 : dof_trans_(&dof_trans)
162 , vdim_(vdim)
163 , ordering_(ordering)
164 {}
165
166 /** @brief Configure the transformation using face orientations for the
167 current element. */
168 /// The face_orientation array can be obtained from Mesh::GetElementFaces.
169 inline void SetFaceOrientations(const Array<int> & Fo)
170 { Fo_ = Fo; }
171
172 /// Return the face orientations for the current element
173 inline const Array<int> & GetFaceOrientations() const { return Fo_; }
174
175 /// Set or change the nested StatelessDofTransformation object
177 {
178 dof_trans_ = &dof_trans;
179 }
181 {
182 dof_trans_ = dof_trans;
183 }
184
185 /// Return the nested StatelessDofTransformation object
187 { return dof_trans_; }
188
189 /// Set or change the vdim and ordering parameter
190 inline void SetVDim(int vdim = 1, int ordering = 0)
191 {
192 vdim_ = vdim;
193 ordering_ = ordering;
194 }
195
196 /// Return the current vdim value
197 inline int GetVDim() const { return vdim_; }
198
199 inline int Size() const { return dof_trans_->Size(); }
200 inline int Height() const { return dof_trans_->Height(); }
201 inline int NumRows() const { return dof_trans_->NumRows(); }
202 inline int Width() const { return dof_trans_->Width(); }
203 inline int NumCols() const { return dof_trans_->NumCols(); }
204 inline bool IsIdentity() const { return dof_trans_->IsIdentity(); }
205
206 /** Transform local DoFs to align with the global DoFs. For example, this
207 transformation can be used to map the local vector computed by
208 FiniteElement::Project() to the transformed vector stored within a
209 GridFunction object. */
210 void TransformPrimal(real_t *v) const;
211 inline void TransformPrimal(Vector &v) const
212 { TransformPrimal(v.GetData()); }
213
214 /// Transform groups of DoFs stored as dense matrices
215 inline void TransformPrimalCols(DenseMatrix &V) const
216 {
217 for (int c=0; c<V.Width(); c++)
218 {
220 }
221 }
222
223 /** Inverse transform local DoFs. Used to transform DoFs from a global vector
224 back to their element-local form. For example, this must be used to
225 transform the vector obtained using GridFunction::GetSubVector before it
226 can be used to compute a local interpolation.
227 */
228 void InvTransformPrimal(real_t *v) const;
229 inline void InvTransformPrimal(Vector &v) const
231
232 /** Transform dual DoFs as computed by a LinearFormIntegrator before summing
233 into a LinearForm object. */
234 void TransformDual(real_t *v) const;
235 inline void TransformDual(Vector &v) const
236 { TransformDual(v.GetData()); }
237
238 /** Inverse Transform dual DoFs */
239 void InvTransformDual(real_t *v) const;
240 inline void InvTransformDual(Vector &v) const
241 { InvTransformDual(v.GetData()); }
242
243 /** Transform a matrix of dual DoFs entries as computed by a
244 BilinearFormIntegrator before summing into a BilinearForm object. */
245 inline void TransformDual(DenseMatrix &V) const
246 {
249 }
250
251 /// Transform rows of a dense matrix containing dual DoFs
252 inline void TransformDualRows(DenseMatrix &V) const
253 {
254 Vector row;
255 for (int r=0; r<V.Height(); r++)
256 {
257 V.GetRow(r, row);
258 TransformDual(row);
259 V.SetRow(r, row);
260 }
261 }
262
263 /// Transform columns of a dense matrix containing dual DoFs
264 inline void TransformDualCols(DenseMatrix &V) const
265 {
266 for (int c=0; c<V.Width(); c++)
267 {
269 }
270 }
271};
272
273/** Transform a matrix of DoFs entries from different finite element spaces as
274 computed by a DiscreteInterpolator before copying into a
275 DiscreteLinearOperator.
276*/
277void TransformPrimal(const DofTransformation *ran_dof_trans,
278 const DofTransformation *dom_dof_trans,
279 DenseMatrix &elmat);
280
281/** Transform a matrix of dual DoFs entries from different finite element spaces
282 as computed by a BilinearFormIntegrator before summing into a
283 MixedBilinearForm object.
284*/
285void TransformDual(const DofTransformation *ran_dof_trans,
286 const DofTransformation *dom_dof_trans,
287 DenseMatrix &elmat);
288
289/** Abstract base class for high-order Nedelec spaces on elements with
290 triangular faces.
291
292 The Nedelec DoFs on the interior of triangular faces come in pairs which
293 share an interpolation point but have different vector directions. These
294 directions depend on the orientation of the face and can therefore differ in
295 neighboring elements. The mapping required to transform these DoFs can be
296 implemented as series of 2x2 linear transformations. The raw data for these
297 linear transformations is stored in the T_data and TInv_data arrays and can
298 be accessed as DenseMatrices using the GetFaceTransform() and
299 GetFaceInverseTransform() methods.
300*/
302{
303private:
304 static const real_t T_data[24];
305 static const real_t TInv_data[24];
306 static const DenseTensor T, TInv;
307
308protected:
309 const int order; // basis function order
310 const int nedofs; // number of DoFs per edge
311 const int nfdofs; // number of DoFs per face
312 const int nedges; // number of edges per element
313 const int nfaces; // number of triangular faces per element
314
315 ND_DofTransformation(int size, int order, int num_edges, int num_tri_faces);
316
317public:
318 // Return the 2x2 transformation operator for the given face orientation
319 static const DenseMatrix & GetFaceTransform(int ori) { return T(ori); }
320
321 // Return the 2x2 inverse transformation operator
322 static const DenseMatrix & GetFaceInverseTransform(int ori)
323 { return TInv(ori); }
324
325 bool IsIdentity() const override { return nfdofs < 2; }
326
327 void TransformPrimal(const Array<int> & Fo, real_t *v) const override;
328 void InvTransformPrimal(const Array<int> & Fo, real_t *v) const override;
329 void TransformDual(const Array<int> & Fo, real_t *v) const override;
330 void InvTransformDual(const Array<int> & Fo, real_t *v) const override;
331};
332
333/// Stateless DoF transformation implementation for the Nedelec basis on
334/// triangles
336{
337public:
341};
342
343/// DoF transformation implementation for the Nedelec basis on tetrahedra
345{
346public:
348 : ND_DofTransformation(order*(order + 2)*(order + 3)/2, order, 6, 4)
349 {}
350};
351
352/// DoF transformation implementation for the Nedelec basis on wedge elements
354{
355public:
357 : ND_DofTransformation(3 * order * ((order + 1) * (order + 2))/2,
358 order, 9, 2)
359 {}
360};
361
362} // namespace mfem
363
364#endif // MFEM_DOFTRANSFORM
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void SetRow(int r, const real_t *row)
void GetColumn(int c, Vector &col) const
void GetRow(int r, Vector &row) const
Rank 3 tensor (array of matrices)
void TransformDual(real_t *v) const
Definition doftrans.cpp:81
void TransformDual(DenseMatrix &V) const
Definition doftrans.hpp:245
void InvTransformPrimal(real_t *v) const
Definition doftrans.cpp:49
DofTransformation(const StatelessDofTransformation &dof_trans, int vdim=1, int ordering=0)
Constructor with a known StatelessDofTransformation.
Definition doftrans.hpp:159
void InvTransformDual(Vector &v) const
Definition doftrans.hpp:240
void TransformPrimalCols(DenseMatrix &V) const
Transform groups of DoFs stored as dense matrices.
Definition doftrans.hpp:215
void TransformPrimal(Vector &v) const
Definition doftrans.hpp:211
const StatelessDofTransformation * dof_trans_
Definition doftrans.hpp:145
void SetFaceOrientations(const Array< int > &Fo)
Configure the transformation using face orientations for the current element.
Definition doftrans.hpp:169
int GetVDim() const
Return the current vdim value.
Definition doftrans.hpp:197
void TransformDual(Vector &v) const
Definition doftrans.hpp:235
void TransformDualRows(DenseMatrix &V) const
Transform rows of a dense matrix containing dual DoFs.
Definition doftrans.hpp:252
const Array< int > & GetFaceOrientations() const
Return the face orientations for the current element.
Definition doftrans.hpp:173
void SetDofTransformation(const StatelessDofTransformation *dof_trans)
Definition doftrans.hpp:180
void SetDofTransformation(const StatelessDofTransformation &dof_trans)
Set or change the nested StatelessDofTransformation object.
Definition doftrans.hpp:176
void InvTransformDual(real_t *v) const
Definition doftrans.cpp:113
void TransformDualCols(DenseMatrix &V) const
Transform columns of a dense matrix containing dual DoFs.
Definition doftrans.hpp:264
void InvTransformPrimal(Vector &v) const
Definition doftrans.hpp:229
const StatelessDofTransformation * GetDofTransformation() const
Return the nested StatelessDofTransformation object.
Definition doftrans.hpp:186
void SetVDim(int vdim=1, int ordering=0)
Set or change the vdim and ordering parameter.
Definition doftrans.hpp:190
bool IsIdentity() const
Definition doftrans.hpp:204
DofTransformation(int vdim=1, int ordering=0)
Default constructor which requires that SetDofTransformation be called before use.
Definition doftrans.hpp:152
void TransformPrimal(real_t *v) const
Definition doftrans.cpp:17
ND_DofTransformation(int size, int order, int num_edges, int num_tri_faces)
Definition doftrans.cpp:203
bool IsIdentity() const override
If the DofTransformation performs no transformation.
Definition doftrans.hpp:325
void InvTransformDual(const Array< int > &Fo, real_t *v) const override
Definition doftrans.cpp:291
void TransformPrimal(const Array< int > &Fo, real_t *v) const override
Definition doftrans.cpp:214
void TransformDual(const Array< int > &Fo, real_t *v) const override
Definition doftrans.cpp:266
static const DenseMatrix & GetFaceInverseTransform(int ori)
Definition doftrans.hpp:322
static const DenseMatrix & GetFaceTransform(int ori)
Definition doftrans.hpp:319
void InvTransformPrimal(const Array< int > &Fo, real_t *v) const override
Definition doftrans.cpp:240
DoF transformation implementation for the Nedelec basis on tetrahedra.
Definition doftrans.hpp:345
DoF transformation implementation for the Nedelec basis on wedge elements.
Definition doftrans.hpp:354
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
virtual void InvTransformPrimal(const Array< int > &face_orientation, real_t *v) const =0
virtual ~StatelessDofTransformation()=default
void InvTransformDual(const Array< int > &face_orientation, Vector &v) const
Definition doftrans.hpp:118
void TransformPrimal(const Array< int > &face_orientation, Vector &v) const
Definition doftrans.hpp:92
void InvTransformPrimal(const Array< int > &face_orientation, Vector &v) const
Definition doftrans.hpp:103
virtual void TransformPrimal(const Array< int > &face_orientation, real_t *v) const =0
virtual void InvTransformDual(const Array< int > &face_orientation, real_t *v) const =0
void TransformDual(const Array< int > &face_orientation, Vector &v) const
Definition doftrans.hpp:111
virtual bool IsIdentity() const =0
If the DofTransformation performs no transformation.
virtual void TransformDual(const Array< int > &face_orientation, real_t *v) const =0
Vector data type.
Definition vector.hpp:80
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
void TransformDual(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition doftrans.cpp:160
void TransformPrimal(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition doftrans.cpp:145
float real_t
Definition config.hpp:43