MFEM  v4.6.0
Finite element discretization library
doftrans.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 
19 namespace 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 {
70 protected:
71  int size_;
72 
74  : size_(size) {}
75 
76 public:
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  /** Transform local DoFs to align with the global DoFs. For example, this
84  transformation can be used to map the local vector computed by
85  FiniteElement::Project() to the transformed vector stored within a
86  GridFunction object. */
87  virtual void TransformPrimal(const Array<int> & face_orientation,
88  double *v) const = 0;
89  inline void TransformPrimal(const Array<int> & face_orientation,
90  Vector &v) const
91  { TransformPrimal(face_orientation, v.GetData()); }
92 
93  /** Inverse transform local DoFs. Used to transform DoFs from a global vector
94  back to their element-local form. For example, this must be used to
95  transform the vector obtained using GridFunction::GetSubVector before it
96  can be used to compute a local interpolation.
97  */
98  virtual void InvTransformPrimal(const Array<int> & face_orientation,
99  double *v) const = 0;
100  inline void InvTransformPrimal(const Array<int> & face_orientation,
101  Vector &v) const
102  { InvTransformPrimal(face_orientation, v.GetData()); }
103 
104  /** Transform dual DoFs as computed by a LinearFormIntegrator before summing
105  into a LinearForm object. */
106  virtual void TransformDual(const Array<int> & face_orientation,
107  double *v) const = 0;
108  inline void TransformDual(const Array<int> & face_orientation,
109  Vector &v) const
110  { TransformDual(face_orientation, v.GetData()); }
111 
112  /** Inverse Transform dual DoFs */
113  virtual void InvTransformDual(const Array<int> & face_orientation,
114  double *v) const = 0;
115  inline void InvTransformDual(const Array<int> & face_orientation,
116  Vector &v) const
117  { InvTransformDual(face_orientation, v.GetData()); }
118 };
119 
120 /** The DofTransformation class is an extension of the
121  StatelessDofTransformation which stores the face orientations used to
122  select the necessary transformations which allows it to offer a collection
123  of convenience methods.
124 
125  DofTransformation objects are provided by the FiniteElementSpace which has
126  access to the mesh and can therefore provide the face orientations. This is
127  convenient when working with GridFunction, LinearForm, or BilinearForm
128  objects or their parallel counterparts.
129 
130  StatelessDofTransformation objects are provided by FiniteElement or
131  FiniteElementCollection objects which do not have access to face
132  orientation information. This can be useful in non-standard contexts such as
133  transferring finite element degrees of freedom between different meshes.
134  For examples of its use see the TransferMap used by the SubMesh class.
135  */
137 {
138 protected:
140 
142  : StatelessDofTransformation(size) {}
143 
144 public:
145 
146  /** @brief Configure the transformation using face orientations for the
147  current element. */
148  /// The face_orientation array can be obtained from Mesh::GetElementFaces.
149  inline void SetFaceOrientations(const Array<int> & face_orientation)
150  { Fo = face_orientation; }
151 
152  inline const Array<int> & GetFaceOrientations() const { return Fo; }
153 
158 
159  /** Transform local DoFs to align with the global DoFs. For example, this
160  transformation can be used to map the local vector computed by
161  FiniteElement::Project() to the transformed vector stored within a
162  GridFunction object. */
163  inline void TransformPrimal(double *v) const
164  { TransformPrimal(Fo, v); }
165  inline void TransformPrimal(Vector &v) const
166  { TransformPrimal(v.GetData()); }
167 
168  /// Transform groups of DoFs stored as dense matrices
169  inline void TransformPrimalCols(DenseMatrix &V) const
170  {
171  for (int c=0; c<V.Width(); c++)
172  {
174  }
175  }
176 
177  /** Inverse transform local DoFs. Used to transform DoFs from a global vector
178  back to their element-local form. For example, this must be used to
179  transform the vector obtained using GridFunction::GetSubVector before it
180  can be used to compute a local interpolation.
181  */
182  inline void InvTransformPrimal(double *v) const
183  { InvTransformPrimal(Fo, v); }
184  inline void InvTransformPrimal(Vector &v) const
185  { InvTransformPrimal(v.GetData()); }
186 
187  /** Transform dual DoFs as computed by a LinearFormIntegrator before summing
188  into a LinearForm object. */
189  inline void TransformDual(double *v) const
190  { TransformDual(Fo, v); }
191  inline void TransformDual(Vector &v) const
192  { TransformDual(v.GetData()); }
193 
194  /** Inverse Transform dual DoFs */
195  inline void InvTransformDual(double *v) const
196  { InvTransformDual(Fo, v); }
197  inline void InvTransformDual(Vector &v) const
198  { InvTransformDual(v.GetData()); }
199 
200  /** Transform a matrix of dual DoFs entries as computed by a
201  BilinearFormIntegrator before summing into a BilinearForm object. */
202  inline void TransformDual(DenseMatrix &V) const
203  {
206  }
207 
208  /// Transform rows of a dense matrix containing dual DoFs
209  inline void TransformDualRows(DenseMatrix &V) const
210  {
211  Vector row;
212  for (int r=0; r<V.Height(); r++)
213  {
214  V.GetRow(r, row);
215  TransformDual(row);
216  V.SetRow(r, row);
217  }
218  }
219 
220  /// Transform columns of a dense matrix containing dual DoFs
221  inline void TransformDualCols(DenseMatrix &V) const
222  {
223  for (int c=0; c<V.Width(); c++)
224  {
225  TransformDual(V.GetColumn(c));
226  }
227  }
228 
229  virtual ~DofTransformation() = default;
230 };
231 
232 /** Transform a matrix of DoFs entries from different finite element spaces as
233  computed by a DiscreteInterpolator before copying into a
234  DiscreteLinearOperator.
235 */
236 void TransformPrimal(const DofTransformation *ran_dof_trans,
237  const DofTransformation *dom_dof_trans,
238  DenseMatrix &elmat);
239 
240 /** Transform a matrix of dual DoFs entries from different finite element spaces
241  as computed by a BilinearFormIntegrator before summing into a
242  MixedBilinearForm object.
243 */
244 void TransformDual(const DofTransformation *ran_dof_trans,
245  const DofTransformation *dom_dof_trans,
246  DenseMatrix &elmat);
247 
248 /** The StatelessVDofTransformation class implements a nested transformation
249  where an arbitrary StatelessDofTransformation is replicated with a
250  vdim >= 1.
251 */
253 {
254 protected:
255  int vdim_;
258 
259 public:
260  /** @brief Default constructor which requires that SetDofTransformation be
261  called before use. */
262  StatelessVDofTransformation(int vdim = 1, int ordering = 0)
264  , vdim_(vdim)
265  , ordering_(ordering)
266  , sdoftrans_(NULL)
267  {}
268 
269  /// Constructor with a known StatelessDofTransformation
271  int vdim = 1,
272  int ordering = 0)
273  : StatelessDofTransformation(vdim * doftrans.Size())
274  , vdim_(vdim)
275  , ordering_(ordering)
276  , sdoftrans_(&doftrans)
277  {}
278 
279  /// Set or change the vdim parameter
280  inline void SetVDim(int vdim)
281  {
282  vdim_ = vdim;
283  if (sdoftrans_)
284  {
285  size_ = vdim_ * sdoftrans_->Size();
286  }
287  }
288 
289  /// Return the current vdim value
290  inline int GetVDim() const { return vdim_; }
291 
292  /// Set or change the nested StatelessDofTransformation object
294  {
295  size_ = vdim_ * doftrans.Size();
296  sdoftrans_ = &doftrans;
297  }
298 
299  /// Return the nested StatelessDofTransformation object
301  { return sdoftrans_; }
302 
307 
308  /** Specializations of these base class methods which account for the vdim
309  and ordering of the full set of DoFs.
310  */
311  void TransformPrimal(const Array<int> & face_ori, double *v) const;
312  void InvTransformPrimal(const Array<int> & face_ori, double *v) const;
313  void TransformDual(const Array<int> & face_ori, double *v) const;
314  void InvTransformDual(const Array<int> & face_ori, double *v) const;
315 };
316 
317 /** The VDofTransformation class implements a nested transformation where an
318  arbitrary DofTransformation is replicated with a vdim >= 1.
319 */
321  public DofTransformation
322 {
323 protected:
325 
326 public:
327  /** @brief Default constructor which requires that SetDofTransformation be
328  called before use. */
329  VDofTransformation(int vdim = 1, int ordering = 0)
331  , StatelessVDofTransformation(vdim, ordering)
332  , DofTransformation(0)
333  , doftrans_(NULL)
334  {}
335 
336  /// Constructor with a known DofTransformation
337  /// @note The face orientations in @a doftrans will be copied into the
338  /// new VDofTransformation object.
339  VDofTransformation(DofTransformation & doftrans, int vdim = 1,
340  int ordering = 0)
341  : StatelessDofTransformation(vdim * doftrans.Size())
342  , StatelessVDofTransformation(doftrans, vdim, ordering)
343  , DofTransformation(vdim * doftrans.Size())
344  , doftrans_(&doftrans)
345  {
347  }
348 
350 
351  /// Set or change the nested DofTransformation object
352  /// @note The face orientations in @a doftrans will be copied into the
353  /// VDofTransformation object.
355  {
356  doftrans_ = &doftrans;
359  }
360 
361  /// Return the nested DofTransformation object
362  inline DofTransformation * GetDofTransformation() const { return doftrans_; }
363 
364  /// Set new face orientations in both the VDofTransformation and the
365  /// DofTransformation contained within (if there is one).
366  inline void SetFaceOrientations(const Array<int> & face_orientation)
367  {
368  DofTransformation::SetFaceOrientations(face_orientation);
369  if (doftrans_) { doftrans_->SetFaceOrientations(face_orientation); }
370  }
371 
376 
377  inline void TransformPrimal(double *v) const
378  { TransformPrimal(Fo, v); }
379  inline void InvTransformPrimal(double *v) const
380  { InvTransformPrimal(Fo, v); }
381  inline void TransformDual(double *v) const
382  { TransformDual(Fo, v); }
383  inline void InvTransformDual(double *v) const
384  { InvTransformDual(Fo, v); }
385 };
386 
387 /** Abstract base class for high-order Nedelec spaces on elements with
388  triangular faces.
389 
390  The Nedelec DoFs on the interior of triangular faces come in pairs which
391  share an interpolation point but have different vector directions. These
392  directions depend on the orientation of the face and can therefore differ in
393  neighboring elements. The mapping required to transform these DoFs can be
394  implemented as series of 2x2 linear transformations. The raw data for these
395  linear transformations is stored in the T_data and TInv_data arrays and can
396  be accessed as DenseMatrices using the GetFaceTransform() and
397  GetFaceInverseTransform() methods.
398 */
400 {
401 private:
402  static const double T_data[24];
403  static const double TInv_data[24];
404  static const DenseTensor T, TInv;
405 
406 protected:
407  const int order; // basis function order
408  const int nedofs; // number of DoFs per edge
409  const int nfdofs; // number of DoFs per face
410  const int nedges; // number of edges per element
411  const int nfaces; // number of triangular faces per element
412 
413  ND_StatelessDofTransformation(int size, int order,
414  int num_edges, int num_tri_faces);
415 
416 public:
417  // Return the 2x2 transformation operator for the given face orientation
418  static const DenseMatrix & GetFaceTransform(int ori) { return T(ori); }
419 
420  // Return the 2x2 inverse transformation operator
421  static const DenseMatrix & GetFaceInverseTransform(int ori)
422  { return TInv(ori); }
423 
424  void TransformPrimal(const Array<int> & face_orientation,
425  double *v) const;
426 
427  void InvTransformPrimal(const Array<int> & face_orientation,
428  double *v) const;
429 
430  void TransformDual(const Array<int> & face_orientation,
431  double *v) const;
432 
433  void InvTransformDual(const Array<int> & face_orientation,
434  double *v) const;
435 };
436 
437 /// Stateless DoF transformation implementation for the Nedelec basis on
438 /// triangles
440 {
441 public:
445  {}
446 };
447 
448 /// DoF transformation implementation for the Nedelec basis on triangles
451 {
452 public:
455  , DofTransformation(order*(order + 2))
457  {}
458 
463 
468 };
469 
470 /// DoF transformation implementation for the Nedelec basis on tetrahedra
472 {
473 public:
477  6, 4)
478  {}
479 };
480 
481 /// DoF transformation implementation for the Nedelec basis on tetrahedra
484 {
485 public:
488  , DofTransformation(order*(order + 2)*(order + 3)/2)
490  {}
491 
496 
501 };
502 
503 /// DoF transformation implementation for the Nedelec basis on wedge elements
505 {
506 public:
508  : StatelessDofTransformation(3 * order * ((order + 1) * (order + 2))/2)
509  , ND_StatelessDofTransformation(3 * order * ((order + 1) * (order + 2))/2,
510  order, 9, 2)
511  {}
512 };
513 
514 /// DoF transformation implementation for the Nedelec basis on wedge elements
517 {
518 public:
520  : StatelessDofTransformation(3 * order * ((order + 1) * (order + 2))/2)
521  , DofTransformation(3 * order * ((order + 1) * (order + 2))/2)
523  {}
524 
529 
534 };
535 
536 } // namespace mfem
537 
538 #endif // MFEM_DOFTRANSFORM
DoF transformation implementation for the Nedelec basis on tetrahedra.
Definition: doftrans.hpp:471
void InvTransformDual(const Array< int > &face_orientation, Vector &v) const
Definition: doftrans.hpp:115
void TransformDual(const Array< int > &face_orientation, Vector &v) const
Definition: doftrans.hpp:108
virtual void TransformDual(const Array< int > &face_orientation, double *v) const =0
void InvTransformPrimal(const Array< int > &face_orientation, Vector &v) const
Definition: doftrans.hpp:100
StatelessVDofTransformation(int vdim=1, int ordering=0)
Default constructor which requires that SetDofTransformation be called before use.
Definition: doftrans.hpp:262
void SetRow(int r, const double *row)
Definition: densemat.cpp:2177
VDofTransformation(int vdim=1, int ordering=0)
Default constructor which requires that SetDofTransformation be called before use.
Definition: doftrans.hpp:329
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
void TransformDual(Vector &v) const
Definition: doftrans.hpp:191
void TransformDual(double *v) const
Definition: doftrans.hpp:189
void TransformPrimalCols(DenseMatrix &V) const
Transform groups of DoFs stored as dense matrices.
Definition: doftrans.hpp:169
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void InvTransformPrimal(double *v) const
Definition: doftrans.hpp:379
void InvTransformDual(const Array< int > &face_ori, double *v) const
Definition: doftrans.cpp:154
DofTransformation(int size)
Definition: doftrans.hpp:141
void InvTransformPrimal(Vector &v) const
Definition: doftrans.hpp:184
DofTransformation * GetDofTransformation() const
Return the nested DofTransformation object.
Definition: doftrans.hpp:362
ND_StatelessDofTransformation(int size, int order, int num_edges, int num_tri_faces)
Definition: doftrans.cpp:212
void InvTransformDual(double *v) const
Definition: doftrans.hpp:383
void TransformPrimal(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition: doftrans.cpp:17
const Array< int > & GetFaceOrientations() const
Definition: doftrans.hpp:152
void InvTransformDual(const Array< int > &face_orientation, double *v) const
Definition: doftrans.cpp:296
void SetDofTransformation(DofTransformation &doftrans)
Definition: doftrans.hpp:354
void InvTransformPrimal(double *v) const
Definition: doftrans.hpp:182
void SetFaceOrientations(const Array< int > &face_orientation)
Definition: doftrans.hpp:366
void SetFaceOrientations(const Array< int > &face_orientation)
Configure the transformation using face orientations for the current element.
Definition: doftrans.hpp:149
void TransformDual(DenseMatrix &V) const
Definition: doftrans.hpp:202
void SetVDim(int vdim)
Set or change the vdim parameter.
Definition: doftrans.hpp:280
virtual void TransformPrimal(const Array< int > &face_orientation, double *v) const =0
virtual void InvTransformPrimal(const Array< int > &face_orientation, double *v) const =0
void GetColumn(int c, Vector &col) const
Definition: densemat.cpp:1330
void TransformDual(double *v) const
Definition: doftrans.hpp:381
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:206
void TransformPrimal(double *v) const
Definition: doftrans.hpp:377
DoF transformation implementation for the Nedelec basis on triangles.
Definition: doftrans.hpp:449
void TransformDual(const Array< int > &face_ori, double *v) const
Definition: doftrans.cpp:124
int GetVDim() const
Return the current vdim value.
Definition: doftrans.hpp:290
void GetRow(int r, Vector &row) const
Definition: densemat.cpp:1314
void TransformPrimal(Vector &v) const
Definition: doftrans.hpp:165
void TransformPrimal(double *v) const
Definition: doftrans.hpp:163
DofTransformation * doftrans_
Definition: doftrans.hpp:324
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
void InvTransformDual(Vector &v) const
Definition: doftrans.hpp:197
void InvTransformPrimal(const Array< int > &face_orientation, double *v) const
Definition: doftrans.cpp:248
static const DenseMatrix & GetFaceTransform(int ori)
Definition: doftrans.hpp:418
static const DenseMatrix & GetFaceInverseTransform(int ori)
Definition: doftrans.hpp:421
StatelessVDofTransformation(StatelessDofTransformation &doftrans, int vdim=1, int ordering=0)
Constructor with a known StatelessDofTransformation.
Definition: doftrans.hpp:270
void InvTransformPrimal(const Array< int > &face_ori, double *v) const
Definition: doftrans.cpp:93
void TransformPrimal(const Array< int > &face_orientation, Vector &v) const
Definition: doftrans.hpp:89
void TransformDual(const Array< int > &face_orientation, double *v) const
Definition: doftrans.cpp:272
void TransformDualCols(DenseMatrix &V) const
Transform columns of a dense matrix containing dual DoFs.
Definition: doftrans.hpp:221
VDofTransformation(DofTransformation &doftrans, int vdim=1, int ordering=0)
Definition: doftrans.hpp:339
void InvTransformDual(double *v) const
Definition: doftrans.hpp:195
DoF transformation implementation for the Nedelec basis on wedge elements.
Definition: doftrans.hpp:504
void TransformDualRows(DenseMatrix &V) const
Transform rows of a dense matrix containing dual DoFs.
Definition: doftrans.hpp:209
Vector data type.
Definition: vector.hpp:58
virtual void InvTransformDual(const Array< int > &face_orientation, double *v) const =0
void TransformDual(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition: doftrans.cpp:40
void TransformPrimal(const Array< int > &face_ori, double *v) const
Definition: doftrans.cpp:63
DoF transformation implementation for the Nedelec basis on tetrahedra.
Definition: doftrans.hpp:482
StatelessDofTransformation * GetDofTransformation() const
Return the nested StatelessDofTransformation object.
Definition: doftrans.hpp:300
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:1096
void TransformPrimal(const Array< int > &face_orientation, double *v) const
Definition: doftrans.cpp:224
void SetDofTransformation(StatelessDofTransformation &doftrans)
Set or change the nested StatelessDofTransformation object.
Definition: doftrans.hpp:293
virtual ~DofTransformation()=default
DoF transformation implementation for the Nedelec basis on wedge elements.
Definition: doftrans.hpp:515
StatelessDofTransformation * sdoftrans_
Definition: doftrans.hpp:257