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