MFEM  v4.5.2
Finite element discretization library
restriction.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_RESTRICTION
13 #define MFEM_RESTRICTION
14 
15 #include "../linalg/operator.hpp"
16 #include "../mesh/mesh.hpp"
17 
18 namespace mfem
19 {
20 
21 class FiniteElementSpace;
22 enum class ElementDofOrdering;
23 
24 /// Abstract base class that defines an interface for element restrictions.
26 {
27 public:
28  /// @brief Add the E-vector degrees of freedom @a x to the L-vector degrees
29  /// of freedom @a y.
30  void AddMultTranspose(const Vector &x, Vector &y,
31  const double a = 1.0) const override = 0;
32 };
33 
34 /// Operator that converts FiniteElementSpace L-vectors to E-vectors.
35 /** Objects of this type are typically created and owned by FiniteElementSpace
36  objects, see FiniteElementSpace::GetElementRestriction(). */
38 {
39 private:
40  /** This number defines the maximum number of elements any dof can belong to
41  for the FillSparseMatrix method. */
42  static const int MaxNbNbr = 16;
43 
44 protected:
46  const int ne;
47  const int vdim;
48  const bool byvdim;
49  const int ndofs;
50  const int dof;
51  const int nedofs;
55 
56  friend class BatchedLORAssembly;
57  friend class BatchedLOR_ADS;
58  friend class BatchedLOR_AMS;
59 
60  /// @name Low-level access to the underlying element-dof mappings
61  ///@{
62  const Array<int> &GatherMap() const { return gather_map; }
63  const Array<int> &Indices() const { return indices; }
64  const Array<int> &Offsets() const { return offsets; }
65  ///@}
66 
67 public:
69  void Mult(const Vector &x, Vector &y) const override;
70  void MultTranspose(const Vector &x, Vector &y) const override;
71  void AddMultTranspose(const Vector &x, Vector &y,
72  const double a = 1.0) const override;
73 
74  /// Compute Mult without applying signs based on DOF orientations.
75  void MultUnsigned(const Vector &x, Vector &y) const;
76  /// Compute MultTranspose without applying signs based on DOF orientations.
77  void MultTransposeUnsigned(const Vector &x, Vector &y) const;
78 
79  /// Compute MultTranspose by setting (rather than adding) element
80  /// contributions; this is a left inverse of the Mult() operation
81  void MultLeftInverse(const Vector &x, Vector &y) const;
82 
83  /// @brief Fills the E-vector y with `boolean` values 0.0 and 1.0 such that each
84  /// each entry of the L-vector is uniquely represented in `y`.
85  /** This means, the sum of the E-vector `y` is equal to the sum of the
86  corresponding L-vector filled with ones. The boolean mask is required to
87  emulate SetSubVector and its transpose on GPUs. This method is running on
88  the host, since the `processed` array requires a large shared memory. */
89  void BooleanMask(Vector& y) const;
90 
91  /// Fill a Sparse Matrix with Element Matrices.
92  void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const;
93 
94  /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
95  given by this ElementRestriction. */
96  int FillI(SparseMatrix &mat) const;
97  /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
98  pattern given by this ElementRestriction, and the values of ea_data. */
99  void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const;
100  /// @private Not part of the public interface (device kernel limitation).
101  ///
102  /// Performs either MultTranspose or AddMultTranspose depending on the
103  /// boolean template parameter @a ADD.
104  template <bool ADD> void TAddMultTranspose(const Vector &x, Vector &y) const;
105 };
106 
107 /// Operator that converts L2 FiniteElementSpace L-vectors to E-vectors.
108 /** Objects of this type are typically created and owned by FiniteElementSpace
109  objects, see FiniteElementSpace::GetElementRestriction(). L-vectors
110  corresponding to grid functions in L2 finite element spaces differ from
111  E-vectors only in the ordering of the degrees of freedom. */
113 {
114  const int ne;
115  const int vdim;
116  const bool byvdim;
117  const int ndof;
118  const int ndofs;
119 public:
121  void Mult(const Vector &x, Vector &y) const override;
122  void MultTranspose(const Vector &x, Vector &y) const override;
123  void AddMultTranspose(const Vector &x, Vector &y,
124  const double a = 1.0) const override;
125  /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
126  given by this ElementRestriction. */
127  void FillI(SparseMatrix &mat) const;
128  /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
129  pattern given by this L2FaceRestriction, and the values of ea_data. */
130  void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const;
131  /// @private Not part of the public interface (device kernel limitation).
132  ///
133  /// Performs either MultTranspose or AddMultTranspose depending on the
134  /// boolean template parameter @a ADD.
135  template <bool ADD> void TAddMultTranspose(const Vector &x, Vector &y) const;
136 };
137 
138 /** An enum type to specify if only e1 value is requested (SingleValued) or both
139  e1 and e2 (DoubleValued). */
141 
142 /** @brief Base class for operators that extracts Face degrees of freedom.
143 
144  In order to compute quantities on the faces of a mesh, it is often useful to
145  extract the degrees of freedom on the faces of the elements. This class
146  provides an interface for such operations.
147 
148  If the FiniteElementSpace is ordered by Ordering::byVDIM, then the expected
149  format for the L-vector is (vdim x ndofs), otherwise if Ordering::byNODES
150  the expected format is (ndofs x vdim), where ndofs is the total number of
151  degrees of freedom.
152  Since FiniteElementSpace can either be continuous or discontinuous, the
153  degrees of freedom on a face can either be single valued or double valued,
154  this is what we refer to as the multiplicity and is represented by the
155  L2FaceValues enum type.
156  The format of the output face E-vector of degrees of freedom is
157  (face_dofs x vdim x multiplicity x nfaces), where face_dofs is the number of
158  degrees of freedom on each face, and nfaces the number of faces of the
159  requested FaceType (see FiniteElementSpace::GetNFbyType).
160 
161  @note Objects of this type are typically created and owned by
162  FiniteElementSpace objects, see FiniteElementSpace::GetFaceRestriction(). */
163 class FaceRestriction : public Operator
164 {
165 public:
167 
168  FaceRestriction(int h, int w): Operator(h, w) { }
169 
170  virtual ~FaceRestriction() { }
171 
172  /** @brief Extract the face degrees of freedom from @a x into @a y.
173 
174  @param[in] x The L-vector of degrees of freedom.
175  @param[out] y The degrees of freedom on the face, corresponding to a face
176  E-vector.
177  */
178  void Mult(const Vector &x, Vector &y) const override = 0;
179 
180  /** @brief Add the face degrees of freedom @a x to the element degrees of
181  freedom @a y.
182 
183  @param[in] x The face degrees of freedom on the face.
184  @param[in,out] y The L-vector of degrees of freedom to which we add the
185  face degrees of freedom.
186  @param[in] a Scalar coefficient for addition.
187  */
188  virtual void AddMultTranspose(const Vector &x, Vector &y,
189  const double a = 1.0) const override = 0;
190 
191  /** @brief Add the face degrees of freedom @a x to the element degrees of
192  freedom @a y. Perform the same computation as AddMultTranspose, but
193  @a x is invalid after calling this method.
194 
195  @param[in,out] x The face degrees of freedom on the face.
196  @param[in,out] y The L-vector of degrees of freedom to which we add the
197  face degrees of freedom.
198 
199  @note This method is an optimization of AddMultTranspose where the @a x
200  Vector is used and modified to avoid memory allocation and memcpy.
201  */
202  virtual void AddMultTransposeInPlace(Vector &x, Vector &y) const
203  {
204  AddMultTranspose(x, y);
205  }
206 
207  /** @brief Set the face degrees of freedom in the element degrees of freedom
208  @a y to the values given in @a x.
209 
210  @param[in] x The face degrees of freedom on the face.
211  @param[in,out] y The L-vector of degrees of freedom to which we add the
212  face degrees of freedom.
213  */
214  void MultTranspose(const Vector &x, Vector &y) const override
215  {
216  y = 0.0;
217  AddMultTranspose(x, y);
218  }
219 };
220 
221 /// Operator that extracts Face degrees of freedom for H1 FiniteElementSpaces.
222 /** Objects of this type are typically created and owned by FiniteElementSpace
223  objects, see FiniteElementSpace::GetFaceRestriction(). */
225 {
226 protected:
228  const int nf; // Number of faces of the requested type
229  const int vdim;
230  const bool byvdim;
231  const int face_dofs; // Number of dofs on each face
232  const int elem_dofs; // Number of dofs in each element
233  const int nfdofs; // Total number of face E-vector dofs
234  const int ndofs; // Total number of dofs
235  Array<int> scatter_indices; // Scattering indices for element 1 on each face
236  Array<int> gather_offsets; // offsets for the gathering indices of each dof
237  Array<int> gather_indices; // gathering indices for each dof
238 
239  /** @brief Construct an H1FaceRestriction.
240 
241  @param[in] fes The FiniteElementSpace on which this operates
242  @param[in] ordering Request a specific element ordering
243  @param[in] type Request internal or boundary faces dofs
244  @param[in] build Request the NCL2FaceRestriction to compute the
245  scatter/gather indices. False should only be used
246  when inheriting from H1FaceRestriction.
247  */
249  const ElementDofOrdering ordering,
250  const FaceType type,
251  bool build);
252 public:
253  /** @brief Construct an H1FaceRestriction.
254 
255  @param[in] fes The FiniteElementSpace on which this operates
256  @param[in] ordering Request a specific element ordering
257  @param[in] type Request internal or boundary faces dofs */
259  const ElementDofOrdering ordering,
260  const FaceType type);
261 
262  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
263  face E-Vector.
264 
265  @param[in] x The L-vector degrees of freedom.
266  @param[out] y The face E-Vector degrees of freedom with the given format:
267  face_dofs x vdim x nf
268  where nf is the number of interior or boundary faces
269  requested by @a type in the constructor.
270  The face_dofs are ordered according to the given
271  ElementDofOrdering. */
272  void Mult(const Vector &x, Vector &y) const override;
273 
275 
276  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
277  L-Vector.
278 
279  @param[in] x The face E-Vector degrees of freedom with the given format:
280  face_dofs x vdim x nf
281  where nf is the number of interior or boundary faces
282  requested by @a type in the constructor.
283  The face_dofs should be ordered according to the given
284  ElementDofOrdering
285  @param[in,out] y The L-vector degrees of freedom.
286  @param[in] a Scalar coefficient for addition. */
287  void AddMultTranspose(const Vector &x, Vector &y,
288  const double a = 1.0) const override;
289 
290 private:
291  /** @brief Compute the scatter indices: L-vector to E-vector, and the offsets
292  for the gathering: E-vector to L-vector.
293 
294  @param[in] ordering Request a specific element ordering.
295  @param[in] type Request internal or boundary faces dofs.
296  */
297  void ComputeScatterIndicesAndOffsets(const ElementDofOrdering ordering,
298  const FaceType type);
299 
300  /** @brief Compute the gather indices: E-vector to L-vector.
301 
302  Note: Requires the gather offsets to be computed.
303 
304  @param[in] ordering Request a specific element ordering.
305  @param[in] type Request internal or boundary faces dofs.
306  */
307  void ComputeGatherIndices(const ElementDofOrdering ordering,
308  const FaceType type);
309 
310 protected:
311  mutable Array<int> face_map; // Used in the computation of GetFaceDofs
312 
313  /** @brief Verify that H1FaceRestriction is build from an H1 FESpace.
314 
315  @param[in] ordering The FESpace element ordering.
316  */
317  void CheckFESpace(const ElementDofOrdering ordering);
318 
319  /** @brief Set the scattering indices of elem1, and increment the offsets for
320  the face described by the @a face.
321 
322  @param[in] face The face information of the current face.
323  @param[in] face_index The interior/boundary face index.
324  @param[in] ordering Request a specific element ordering.
325  */
327  const int face_index,
328  const ElementDofOrdering ordering);
329 
330  /** @brief Set the gathering indices of elem1 for the interior face described
331  by the @a face.
332 
333  @param[in] face The face information of the current face.
334  @param[in] face_index The interior/boundary face index.
335  @param[in] ordering Request a specific element ordering.
336  */
338  const int face_index,
339  const ElementDofOrdering ordering);
340 };
341 
342 /// Operator that extracts Face degrees of freedom for L2 spaces.
343 /** Objects of this type are typically created and owned by FiniteElementSpace
344  objects, see FiniteElementSpace::GetFaceRestriction(). */
346 {
347 protected:
349  const int nf; // Number of faces of the requested type
350  const int ne; // Number of elements
351  const int vdim; // vdim
352  const bool byvdim;
353  const int face_dofs; // Number of dofs on each face
354  const int elem_dofs; // Number of dofs in each element
355  const int nfdofs; // Total number of dofs on the faces
356  const int ndofs; // Total number of dofs
357  const FaceType type;
359  Array<int> scatter_indices1; // Scattering indices for element 1 on each face
360  Array<int> scatter_indices2; // Scattering indices for element 2 on each face
361  Array<int> gather_offsets; // offsets for the gathering indices of each dof
362  Array<int> gather_indices; // gathering indices for each dof
363 
364  /** @brief Constructs an L2FaceRestriction.
365 
366  @param[in] fes The FiniteElementSpace on which this operates
367  @param[in] ordering Request a specific ordering
368  @param[in] type Request internal or boundary faces dofs
369  @param[in] m Request the face dofs for elem1, or both elem1 and
370  elem2
371  @param[in] build Request the NCL2FaceRestriction to compute the
372  scatter/gather indices. False should only be used
373  when inheriting from L2FaceRestriction.
374  */
376  const ElementDofOrdering ordering,
377  const FaceType type,
378  const L2FaceValues m,
379  bool build);
380 
381 public:
382  /** @brief Constructs an L2FaceRestriction.
383 
384  @param[in] fes The FiniteElementSpace on which this operates
385  @param[in] ordering Request a specific ordering
386  @param[in] type Request internal or boundary faces dofs
387  @param[in] m Request the face dofs for elem1, or both elem1 and
388  elem2 */
390  const ElementDofOrdering ordering,
391  const FaceType type,
393 
394  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
395  face E-Vector.
396 
397  @param[in] x The L-vector degrees of freedom.
398  @param[out] y The face E-Vector degrees of freedom with the given format:
399  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf)
400  if L2FacesValues::SingleValued (face_dofs x vdim x nf)
401  where nf is the number of interior or boundary faces
402  requested by @a type in the constructor.
403  The face_dofs are ordered according to the given
404  ElementDofOrdering. */
405  void Mult(const Vector &x, Vector &y) const override;
406 
408 
409  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
410  L-Vector.
411 
412  @param[in] x The face E-Vector degrees of freedom with the given format:
413  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf)
414  if L2FacesValues::SingleValued (face_dofs x vdim x nf)
415  where nf is the number of interior or boundary faces
416  requested by @a type in the constructor.
417  The face_dofs should be ordered according to the given
418  ElementDofOrdering
419  @param[in,out] y The L-vector degrees of freedom.
420  @param[in] a Scalar coefficient for addition. */
421  void AddMultTranspose(const Vector &x, Vector &y,
422  const double a = 1.0) const override;
423 
424  /** @brief Fill the I array of SparseMatrix corresponding to the sparsity
425  pattern given by this L2FaceRestriction.
426 
427  @param[in,out] mat The sparse matrix for which we want to initialize the
428  row offsets.
429  @param[in] keep_nbr_block When set to true the SparseMatrix will
430  include the rows (in addition to the columns)
431  corresponding to face-neighbor dofs. The
432  default behavior is to disregard those rows. */
433  virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block = false) const;
434 
435  /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
436  the sparsity pattern given by this L2FaceRestriction, and the values of
437  fea_data.
438 
439  @param[in] fea_data The dense matrices representing the local operators
440  on each face. The format is:
441  face_dofs x face_dofs x 2 x nf
442  On each face the first local matrix corresponds to
443  the contribution of elem1 on elem2, and the second to
444  the contribution of elem2 on elem1.
445  @param[in,out] mat The sparse matrix that is getting filled.
446  @param[in] keep_nbr_block When set to true the SparseMatrix will
447  include the rows (in addition to the columns)
448  corresponding to face-neighbor dofs. The
449  default behavior is to disregard those rows. */
450  virtual void FillJAndData(const Vector &fea_data,
451  SparseMatrix &mat,
452  const bool keep_nbr_block = false) const;
453 
454  /** @brief This methods adds the DG face matrices to the element matrices.
455 
456  @param[in] fea_data The dense matrices representing the local operators
457  on each face. The format is:
458  face_dofs x face_dofs x 2 x nf
459  On each face the first and second local matrices
460  correspond to the contributions of elem1 and elem2 on
461  themselves respectively.
462  @param[in,out] ea_data The dense matrices representing the element local
463  contributions for each element to which will be
464  added the face contributions.
465  The format is: dofs x dofs x ne, where dofs is the
466  number of dofs per element and ne the number of
467  elements. */
468  virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data,
469  Vector &ea_data) const;
470 
471 private:
472  /** @brief Compute the scatter indices: L-vector to E-vector, and the offsets
473  for the gathering: E-vector to L-vector.
474 
475  @param[in] ordering Request a specific element ordering.
476  @param[in] type Request internal or boundary faces dofs.
477  */
478  void ComputeScatterIndicesAndOffsets(const ElementDofOrdering ordering,
479  const FaceType type);
480 
481  /** @brief Compute the gather indices: E-vector to L-vector.
482 
483  Note: Requires the gather offsets to be computed.
484 
485  @param[in] ordering Request a specific element ordering.
486  @param[in] type Request internal or boundary faces dofs.
487  */
488  void ComputeGatherIndices(const ElementDofOrdering ordering,
489  const FaceType type);
490 
491 protected:
492  mutable Array<int> face_map; // Used in the computation of GetFaceDofs
493 
494  /** @brief Verify that L2FaceRestriction is build from an L2 FESpace.
495 
496  @param[in] ordering The FESpace element ordering.
497  */
498  void CheckFESpace(const ElementDofOrdering ordering);
499 
500  /** @brief Set the scattering indices of elem1, and increment the offsets for
501  the face described by the @a face. The ordering of the face dofs of elem1
502  is lexicographic relative to elem1.
503 
504  @param[in] face The face information of the current face.
505  @param[in] face_index The interior/boundary face index.
506  */
508  const int face_index);
509 
510  /** @brief Permute and set the scattering indices of elem2, and increment the
511  offsets for the face described by the @a face. The permutation orders the
512  dofs of elem2 lexicographically as the ones of elem1.
513 
514  @param[in] face The face information of the current face.
515  @param[in] face_index The interior/boundary face index.
516  */
518  const int face_index);
519 
520  /** @brief Permute and set the scattering indices of elem2 for the shared
521  face described by the @a face. The permutation orders the dofs of elem2 as
522  the ones of elem1.
523 
524  @param[in] face The face information of the current face.
525  @param[in] face_index The interior/boundary face index.
526  */
528  const Mesh::FaceInformation &face,
529  const int face_index);
530 
531  /** @brief Set the scattering indices of elem2 for the boundary face
532  described by the @a face.
533 
534  @param[in] face The face information of the current face.
535  @param[in] face_index The interior/boundary face index.
536  */
538  const int face_index);
539 
540  /** @brief Set the gathering indices of elem1 for the interior face described
541  by the @a face.
542 
543  Note: This function modifies the offsets.
544 
545  @param[in] face The face information of the current face.
546  @param[in] face_index The interior/boundary face index.
547  */
549  const int face_index);
550 
551  /** @brief Permute and set the gathering indices of elem2 for the interior
552  face described by the @a face. The permutation orders the dofs of elem2 as
553  the ones of elem1.
554 
555  Note: This function modifies the offsets.
556 
557  @param[in] face The face information of the current face.
558  @param[in] face_index The interior/boundary face index.
559  */
561  const int face_index);
562 
563 public:
564  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
565  face E-Vector. Should only be used with conforming faces and when:
566  m == L2FacesValues::SingleValued
567 
568  @param[in] x The L-vector degrees of freedom.
569  @param[out] y The face E-Vector degrees of freedom with the given format:
570  face_dofs x vdim x nf
571  where nf is the number of interior or boundary faces
572  requested by @a type in the constructor.
573  The face_dofs are ordered according to the given
574  ElementDofOrdering. */
575  void SingleValuedConformingMult(const Vector& x, Vector& y) const;
576 
577  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
578  face E-Vector. Should only be used with conforming faces and when:
579  m == L2FacesValues::DoubleValued
580 
581  @param[in] x The L-vector degrees of freedom.
582  @param[out] y The face E-Vector degrees of freedom with the given format:
583  face_dofs x vdim x 2 x nf
584  where nf is the number of interior or boundary faces
585  requested by @a type in the constructor.
586  The face_dofs are ordered according to the given
587  ElementDofOrdering. */
588  virtual void DoubleValuedConformingMult(const Vector& x, Vector& y) const;
589 
590  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
591  L-Vector. Should only be used with conforming faces and when:
592  m == L2FacesValues::SingleValued
593 
594  @param[in] x The face E-Vector degrees of freedom with the given format:
595  face_dofs x vdim x nf
596  where nf is the number of interior or boundary faces
597  requested by @a type in the constructor.
598  The face_dofs should be ordered according to the given
599  ElementDofOrdering
600  @param[in,out] y The L-vector degrees of freedom. */
601  void SingleValuedConformingAddMultTranspose(const Vector& x, Vector& y) const;
602 
603  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
604  L-Vector. Should only be used with conforming faces and when:
605  m == L2FacesValues::DoubleValued
606 
607  @param[in] x The face E-Vector degrees of freedom with the given format:
608  face_dofs x vdim x 2 x nf
609  where nf is the number of interior or boundary faces
610  requested by @a type in the constructor.
611  The face_dofs should be ordered according to the given
612  ElementDofOrdering
613  @param[in,out] y The L-vector degrees of freedom. */
614  void DoubleValuedConformingAddMultTranspose(const Vector& x, Vector& y) const;
615 };
616 
617 /** This struct stores which side is the master nonconforming side and the
618  index of the interpolator, see InterpolationManager class below. */
620 {
621  uint32_t is_non_conforming : 1;
622  uint32_t master_side : 1;
623  uint32_t index : 30;
624 
625  // default constructor, create a conforming face with index 0.
626  InterpConfig() = default;
627 
628  // Non-conforming face
629  InterpConfig(int master_side, int nc_index)
631  { }
632 
633  InterpConfig(const InterpConfig&) = default;
634 
635  InterpConfig &operator=(const InterpConfig &rhs) = default;
636 };
637 
638 /** This struct stores which side is the master nonconforming side and the
639  index of the interpolator, see InterpolationManager class below. */
641 {
643  uint32_t is_non_conforming : 1;
644  uint32_t master_side : 1;
645  uint32_t index : 30;
646 
647  // default constructor.
648  NCInterpConfig() = default;
649 
650  // Non-conforming face
651  NCInterpConfig(int face_index, int master_side, int nc_index)
655  index(nc_index)
656  { }
657 
658  // Non-conforming face
662  master_side(config.master_side),
663  index(config.index)
664  { }
665 
666  NCInterpConfig(const NCInterpConfig&) = default;
667 
668  NCInterpConfig &operator=(const NCInterpConfig &rhs) = default;
669 };
670 
671 /** @brief This class manages the storage and computation of the interpolations
672  from master (coarse) face to slave (fine) face.
673 */
675 {
676 protected:
679  Array<InterpConfig> interp_config; // interpolator index for each face
680  Array<NCInterpConfig> nc_interp_config; // interpolator index for each ncface
681  Vector interpolators; // face_dofs x face_dofs x num_interpolators
682  int nc_cpt; // Counter for interpolators, and used as index.
683 
684  /** The interpolators are associated to a key of containing the address of
685  PointMatrix and a local face identifier. */
686  using Key = std::pair<const DenseMatrix*,int>;
687  /// The temporary map used to store the different interpolators.
688  using Map = std::map<Key, std::pair<int,const DenseMatrix*>>;
689  Map interp_map; // The temporary map that stores the interpolators.
690 
691 public:
692  InterpolationManager() = delete;
693 
694  /** @brief main constructor.
695 
696  @param[in] fes The FiniteElementSpace on which this operates
697  @param[in] ordering Request a specific element ordering.
698  @param[in] type Request internal or boundary faces dofs
699  */
702  FaceType type);
703 
704  /** @brief Register the face with @a face and index @a face_index as a
705  conforming face for the interpolation of the degrees of freedom.
706 
707  @param[in] face The face information of the current face.
708  @param[in] face_index The interior/boundary face index.
709  */
711  int face_index);
712 
713  /** @brief Register the face with @a face and index @a face_index as a
714  conforming face for the interpolation of the degrees of freedom.
715 
716  @param[in] face The face information of the current face.
717  @param[in] face_index The interior/boundary face index.
718  */
720  int face_index);
721 
722  /** @brief Transform the interpolation matrix map into a contiguous memory
723  structure. */
725 
727 
728  /// @brief Return the total number of interpolators.
730  {
731  return nc_cpt;
732  }
733 
734  /** @brief Return an mfem::Vector containing the interpolators in the
735  following format: face_dofs x face_dofs x num_interpolators. */
736  const Vector& GetInterpolators() const
737  {
738  return interpolators;
739  }
740 
741  /** @brief Return an array containing the interpolation configuration for
742  each face registered with RegisterFaceConformingInterpolation and
743  RegisterFaceCoarseToFineInterpolation. */
745  {
746  return interp_config;
747  }
748 
749  /** @brief Return an array containing the interpolation configuration for
750  each face registered with RegisterFaceConformingInterpolation and
751  RegisterFaceCoarseToFineInterpolation. */
753  {
754  return nc_interp_config;
755  }
756 
757 private:
758  /** @brief Returns the interpolation operator from a master (coarse) face to
759  a slave (fine) face.
760 
761  @param[in] face The face information of the current face.
762  @param[in] ptMat The PointMatrix describing the position and orientation
763  of the fine face in the coarse face. This PointMatrix is
764  usually obtained from the mesh through the method
765  GetNCFacesPtMat.
766  @param[in] ordering Request a specific element ordering.
767  @return The dense matrix corresponding to the interpolation of the face
768  degrees of freedom of the master (coarse) face to the slave
769  (fine) face. */
770  const DenseMatrix* GetCoarseToFineInterpolation(
771  const Mesh::FaceInformation &face,
772  const DenseMatrix* ptMat);
773 };
774 
775 /** @brief Operator that extracts face degrees of freedom for L2 nonconforming
776  spaces.
777 
778  In order to support face restrictions on nonconforming meshes, this
779  operator interpolates master (coarse) face degrees of freedom onto the
780  slave (fine) face. This allows face integrators to treat nonconforming
781  faces just as regular conforming faces. */
783 {
784 protected:
786  mutable Vector x_interp;
787 
788  /** @brief Constructs an NCL2FaceRestriction, this is a specialization of a
789  L2FaceRestriction for nonconforming meshes.
790 
791  @param[in] fes The FiniteElementSpace on which this operates
792  @param[in] ordering Request a specific ordering
793  @param[in] type Request internal or boundary faces dofs
794  @param[in] m Request the face dofs for elem1, or both elem1 and
795  elem2
796  @param[in] build Request the NCL2FaceRestriction to compute the
797  scatter/gather indices. False should only be used
798  when inheriting from NCL2FaceRestriction.
799  */
801  const ElementDofOrdering ordering,
802  const FaceType type,
803  const L2FaceValues m,
804  bool build);
805 public:
806  /** @brief Constructs an NCL2FaceRestriction, this is a specialization of a
807  L2FaceRestriction for nonconforming meshes.
808 
809  @param[in] fes The FiniteElementSpace on which this operates
810  @param[in] ordering Request a specific ordering
811  @param[in] type Request internal or boundary faces dofs
812  @param[in] m Request the face dofs for elem1, or both elem1 and
813  elem2
814  */
816  const ElementDofOrdering ordering,
817  const FaceType type,
819 
820  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
821  face E-Vector.
822 
823  @param[in] x The L-vector degrees of freedom.
824  @param[out] y The face E-Vector degrees of freedom with the given format:
825  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
826  if L2FacesValues::SingleValued (face_dofs x vdim x nf),
827  where nf is the number of interior or boundary faces
828  requested by @a type in the constructor.
829  The face_dofs are ordered according to the given
830  ElementDofOrdering. */
831  void Mult(const Vector &x, Vector &y) const override;
832 
833  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
834  L-Vector.
835 
836  @param[in] x The face E-Vector degrees of freedom with the given format:
837  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
838  if L2FacesValues::SingleValued (face_dofs x vdim x nf),
839  where nf is the number of interior or boundary faces
840  requested by @a type in the constructor.
841  The face_dofs should be ordered according to the given
842  ElementDofOrdering
843  @param[in,out] y The L-vector degrees of freedom.
844  @param[in] a Scalar coefficient for addition. */
845  void AddMultTranspose(const Vector &x, Vector &y,
846  const double a = 1.0) const override;
847 
848  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
849  L-Vector.
850 
851  @param[in,out] x The face E-Vector degrees of freedom with the given format:
852  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
853  if L2FacesValues::SingleValued (face_dofs x vdim x nf),
854  where nf is the number of interior or boundary faces
855  requested by @a type in the constructor.
856  The face_dofs should be ordered according to the given
857  ElementDofOrdering
858  @param[in,out] y The L-vector degrees of freedom.
859 
860  @note This method is an optimization of AddMultTranspose where the @a x
861  Vector is used and modified to avoid memory allocation and memcpy. */
862  void AddMultTransposeInPlace(Vector &x, Vector &y) const override;
863 
864  /** @brief Fill the I array of SparseMatrix corresponding to the sparsity
865  pattern given by this NCL2FaceRestriction.
866 
867  @param[in,out] mat The sparse matrix for which we want to initialize the
868  row offsets.
869  @param[in] keep_nbr_block When set to true the SparseMatrix will
870  include the rows (in addition to the columns)
871  corresponding to face-neighbor dofs. The
872  default behavior is to disregard those rows.
873 
874  @warning This method is not implemented yet. */
875  void FillI(SparseMatrix &mat,
876  const bool keep_nbr_block = false) const override;
877 
878  /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
879  the sparsity pattern given by this NCL2FaceRestriction, and the values of
880  ea_data.
881 
882  @param[in] fea_data The dense matrices representing the local operators
883  on each face. The format is:
884  face_dofs x face_dofs x 2 x nf.
885  On each face the first local matrix corresponds to
886  the contribution of elem1 on elem2, and the second to
887  the contribution of elem2 on elem1.
888  @param[in,out] mat The sparse matrix that is getting filled.
889  @param[in] keep_nbr_block When set to true the SparseMatrix will
890  include the rows (in addition to the columns)
891  corresponding to face-neighbor dofs. The
892  default behavior is to disregard those rows.
893 
894  @warning This method is not implemented yet. */
895  void FillJAndData(const Vector &fea_data,
896  SparseMatrix &mat,
897  const bool keep_nbr_block = false) const override;
898 
899  /** @brief This methods adds the DG face matrices to the element matrices.
900 
901  @param[in] fea_data The dense matrices representing the local operators
902  on each face. The format is:
903  face_dofs x face_dofs x 2 x nf.
904  On each face the first and second local matrices
905  correspond to the contributions of elem1 and elem2 on
906  themselves respectively.
907  @param[in,out] ea_data The dense matrices representing the element local
908  contributions for each element to which will be
909  added the face contributions.
910  The format is: dofs x dofs x ne, where dofs is the
911  number of dofs per element and ne the number of
912  elements.
913 
914  @warning This method is not implemented yet. */
915  void AddFaceMatricesToElementMatrices(const Vector &fea_data,
916  Vector &ea_data) const override;
917 
918 private:
919  /** @brief Compute the scatter indices: L-vector to E-vector, the offsets
920  for the gathering: E-vector to L-vector, and the interpolators from
921  coarse to fine face for master non-comforming faces.
922 
923  @param[in] ordering Request a specific element ordering.
924  @param[in] type Request internal or boundary faces dofs.
925  */
926  void ComputeScatterIndicesAndOffsets(const ElementDofOrdering ordering,
927  const FaceType type);
928 
929  /** @brief Compute the gather indices: E-vector to L-vector.
930 
931  Note: Requires the gather offsets to be computed.
932 
933  @param[in] ordering Request a specific element ordering.
934  @param[in] type Request internal or boundary faces dofs.
935  */
936  void ComputeGatherIndices(const ElementDofOrdering ordering,
937  const FaceType type);
938 
939 public:
940  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
941  face E-Vector. Should only be used with nonconforming faces and when:
942  L2FaceValues m == L2FaceValues::DoubleValued
943 
944  @param[in] x The L-vector degrees of freedom.
945  @param[out] y The face E-Vector degrees of freedom with the given format:
946  (face_dofs x vdim x 2 x nf),
947  where nf is the number of interior or boundary faces
948  requested by @a type in the constructor.
949  The face_dofs are ordered according to the given
950  ElementDofOrdering. */
951  virtual void DoubleValuedNonconformingMult(const Vector& x, Vector& y) const;
952 
953  /** @brief Apply a change of basis from coarse element basis to fine element
954  basis for the coarse face dofs.
955 
956  @param[in,out] x The dofs vector that needs coarse dofs to be express in
957  term of the fine basis.
958  */
960 
961  /** @brief Apply a change of basis from fine element basis to coarse element
962  basis for the coarse face dofs. Should only be used when:
963  L2FaceValues m == L2FaceValues::SingleValued
964 
965  @param[in] x The dofs vector that needs coarse dofs to be express in term
966  of the coarse basis, the result is stored in x_interp.
967  */
969 
970  /** @brief Apply a change of basis from fine element basis to coarse element
971  basis for the coarse face dofs. Should only be used when:
972  L2FaceValues m == L2FaceValues::SingleValued
973 
974  @param[in,out] x The dofs vector that needs coarse dofs to be express in
975  term of the coarse basis, the result is stored in x.
976  */
978 
979  /** @brief Apply a change of basis from fine element basis to coarse element
980  basis for the coarse face dofs. Should only be used when:
981  L2FaceValues m == L2FaceValues::DoubleValued
982 
983  @param[in] x The dofs vector that needs coarse dofs to be express in term
984  of the coarse basis, the result is stored in x_interp.
985  */
987 
988  /** @brief Apply a change of basis from fine element basis to coarse element
989  basis for the coarse face dofs. Should only be used when:
990  L2FaceValues m == L2FaceValues::DoubleValued
991 
992  @param[in,out] x The dofs vector that needs coarse dofs to be express in
993  term of the coarse basis, the result is stored in
994  x.
995  */
997 };
998 
999 /** @brief Return the face map that extracts the degrees of freedom for the
1000  requested local face of a quad or hex, returned in Lexicographic order.
1001 
1002  @param[in] dim The dimension of the space
1003  @param[in] face_id The local face identifier
1004  @param[in] dof1d The 1D number of degrees of freedom for each dimension
1005  @param[out] face_map The map that maps each face dof to an element dof
1006 */
1007 void GetFaceDofs(const int dim, const int face_id,
1008  const int dof1d, Array<int> &face_map);
1009 
1010 /** @brief Convert a dof face index from Native ordering to lexicographic
1011  ordering for quads and hexes.
1012 
1013  @param[in] dim The dimension of the element, 2 for quad, 3 for hex
1014  @param[in] face_id The local face identifier
1015  @param[in] size1d The 1D number of degrees of freedom for each dimension
1016  @param[in] index The native index on the face
1017  @return The lexicographic index on the face
1018 */
1019 int ToLexOrdering(const int dim, const int face_id, const int size1d,
1020  const int index);
1021 
1022 /** @brief Compute the dof face index of elem2 corresponding to the given dof
1023  face index.
1024 
1025  @param[in] dim The dimension of the element, 2 for quad, 3 for hex
1026  @param[in] face_id1 The local face identifier of elem1
1027  @param[in] face_id2 The local face identifier of elem2
1028  @param[in] orientation The orientation of elem2 relative to elem1 on the
1029  face
1030  @param[in] size1d The 1D number of degrees of freedom for each dimension
1031  @param[in] index The dof index on elem1
1032  @return The dof index on elem2 facing the dof on elem1
1033 */
1034 int PermuteFaceL2(const int dim, const int face_id1,
1035  const int face_id2, const int orientation,
1036  const int size1d, const int index);
1037 
1038 }
1039 
1040 #endif // MFEM_RESTRICTION
void FillI(SparseMatrix &mat) const
Array< NCInterpConfig > nc_interp_config
const Array< NCInterpConfig > & GetNCFaceInterpConfig() const
Return an array containing the interpolation configuration for each face registered with RegisterFace...
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
H1FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering ordering, const FaceType type, bool build)
Construct an H1FaceRestriction.
void MultTranspose(const Vector &x, Vector &y) const override
Set the face degrees of freedom in the element degrees of freedom y to the values given in x...
void Mult(const Vector &x, Vector &y) const override=0
Extract the face degrees of freedom from x into y.
void SingleValuedNonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
void PermuteAndSetSharedFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2 for the shared face described by the face...
Operator that extracts Face degrees of freedom for H1 FiniteElementSpaces.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
NCInterpConfig(int face_index, int master_side, int nc_index)
InterpConfig & operator=(const InterpConfig &rhs)=default
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
Array< int > gather_indices
Efficient batched assembly of LOR discretizations on device.
Definition: lor_batched.hpp:33
void SetBoundaryDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem2 for the boundary face described by the face.
std::pair< const DenseMatrix *, int > Key
void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const override
This methods adds the DG face matrices to the element matrices.
virtual void DoubleValuedNonconformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with no...
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
void DoubleValuedNonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
NCInterpConfig(int face_index, InterpConfig &config)
virtual void FillJAndData(const Vector &fea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the J and Data arrays of the SparseMatrix corresponding to the sparsity pattern given by this L2...
void SetFaceDofsGatherIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the gathering indices of elem1 for the interior face described by the face.
FaceRestriction(int h, int w)
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override=0
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const override
Fill the I array of SparseMatrix corresponding to the sparsity pattern given by this NCL2FaceRestrict...
const L2FaceValues m
Array< int > gather_offsets
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
Array< int > gather_indices
void PermuteAndSetFaceDofsGatherIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the gathering indices of elem2 for the interior face described by the face...
L2ElementRestriction(const FiniteElementSpace &)
void DoubleValuedNonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs...
Array< InterpConfig > interp_config
InterpConfig(int master_side, int nc_index)
const Array< int > & GatherMap() const
Definition: restriction.hpp:62
void SingleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
void FillJAndData(const Vector &fea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const override
Fill the J and Data arrays of the SparseMatrix corresponding to the sparsity pattern given by this NC...
virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
FaceType
Definition: mesh.hpp:45
void RegisterFaceConformingInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
const FiniteElementSpace & fes
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 ...
virtual void AddMultTransposeInPlace(Vector &x, Vector &y) const
Add the face degrees of freedom x to the element degrees of freedom y. Perform the same computation a...
void PermuteAndSetFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2, and increment the offsets for the face described by ...
Data type sparse matrix.
Definition: sparsemat.hpp:50
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override=0
Add the face degrees of freedom x to the element degrees of freedom y.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
Operator that extracts Face degrees of freedom for L2 spaces.
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void CheckFESpace(const ElementDofOrdering ordering)
Verify that H1FaceRestriction is build from an H1 FESpace.
const FiniteElementSpace & fes
void LinearizeInterpolatorMapIntoVector()
Transform the interpolation matrix map into a contiguous memory structure.
const Array< int > & Offsets() const
Definition: restriction.hpp:64
void MultLeftInverse(const Vector &x, Vector &y) const
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the I array of SparseMatrix corresponding to the sparsity pattern given by this L2FaceRestrictio...
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
Definition: restriction.cpp:27
void GetFaceDofs(const int dim, const int face_id, const int dof1d, Array< int > &face_map)
Return the face map that extracts the degrees of freedom for the requested local face of a quad or he...
This class manages the storage and computation of the interpolations from master (coarse) face to sla...
Abstract base class that defines an interface for element restrictions.
Definition: restriction.hpp:25
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
This structure is used as a human readable output format that decipheres the information contained in...
Definition: mesh.hpp:1417
const ElementDofOrdering ordering
int ToLexOrdering(const int dim, const int face_id, const int size1d, const int index)
Convert a dof face index from Native ordering to lexicographic ordering for quads and hexes...
const FiniteElementSpace & fes
Definition: restriction.hpp:45
void SetFaceDofsGatherIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering ordering)
Set the gathering indices of elem1 for the interior face described by the face.
Array< int > gather_offsets
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 ...
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Definition: restriction.hpp:37
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
virtual void DoubleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
const Array< InterpConfig > & GetFaceInterpConfig() const
Return an array containing the interpolation configuration for each face registered with RegisterFace...
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
NCL2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an NCL2FaceRestriction, this is a specialization of a L2FaceRestriction for nonconforming ...
void MultUnsigned(const Vector &x, Vector &y) const
Compute Mult without applying signs based on DOF orientations.
std::map< Key, std::pair< int, const DenseMatrix * > > Map
The temporary map used to store the different interpolators.
void SingleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
void CheckFESpace(const ElementDofOrdering ordering)
Verify that L2FaceRestriction is build from an L2 FESpace.
Array< int > scatter_indices1
double a
Definition: lissajous.cpp:41
L2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an L2FaceRestriction.
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
int dim
Definition: ex24.cpp:53
int GetNumInterpolators() const
Return the total number of interpolators.
int FillI(SparseMatrix &mat) const
int index(int i, int j, int nx, int ny)
Definition: life.cpp:235
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
const FiniteElementSpace & fes
void DoubleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
void SetFaceDofsScatterIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering ordering)
Set the scattering indices of elem1, and increment the offsets for the face described by the face...
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void SingleValuedNonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
Vector data type.
Definition: vector.hpp:60
int PermuteFaceL2(const int dim, const int face_id1, const int face_id2, const int orientation, const int size1d, const int index)
Compute the dof face index of elem2 corresponding to the given dof face index.
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void SetFaceDofsScatterIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem1, and increment the offsets for the face described by the face...
const Array< int > & Indices() const
Definition: restriction.hpp:63
void RegisterFaceCoarseToFineInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
Base class for operators that extracts Face degrees of freedom.
Abstract operator.
Definition: operator.hpp:24
uint32_t is_non_conforming
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
const Vector & GetInterpolators() const
Return an mfem::Vector containing the interpolators in the following format: face_dofs x face_dofs x ...
Operator that converts L2 FiniteElementSpace L-vectors to E-vectors.
Array< int > scatter_indices
InterpolationManager interpolations
NCInterpConfig & operator=(const NCInterpConfig &rhs)=default
InterpConfig()=default
Array< int > scatter_indices2
void DoubleValuedNonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
void BooleanMask(Vector &y) const
Fills the E-vector y with boolean values 0.0 and 1.0 such that each each entry of the L-vector is uni...
NCInterpConfig()=default