MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
restriction.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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
16#include "../mesh/mesh.hpp"
18
19namespace mfem
20{
21
22class FiniteElementSpace;
23enum class ElementDofOrdering;
24
25class FaceQuadratureSpace;
26
27/// Abstract base class that defines an interface for element restrictions.
29{
30public:
31 /// @brief Add the E-vector degrees of freedom @a x to the L-vector degrees
32 /// of freedom @a y.
33 void AddMultTranspose(const Vector &x, Vector &y,
34 const real_t a = 1.0) const override = 0;
35};
36
37/// Operator that converts FiniteElementSpace L-vectors to E-vectors.
38/** Objects of this type are typically created and owned by FiniteElementSpace
39 objects, see FiniteElementSpace::GetElementRestriction(). */
41{
42protected:
44 const int ne;
45 const int vdim;
46 const bool byvdim;
47 const int ndofs;
48 const int dof;
49 const int nedofs;
53
54public:
56 void Mult(const Vector &x, Vector &y) const override;
57 void MultTranspose(const Vector &x, Vector &y) const override;
58 void AddMultTranspose(const Vector &x, Vector &y,
59 const real_t a = 1.0) const override;
60
61 /// Compute Mult without applying signs based on DOF orientations.
62 void AbsMult(const Vector &x, Vector &y) const override;
63
64 /// Compute MultTranspose without applying signs based on DOF orientations.
65 void AbsMultTranspose(const Vector &x, Vector &y) const override;
66
67 /// @deprecated Use AbsMult() instead.
68 MFEM_DEPRECATED void MultUnsigned(const Vector &x, Vector &y) const
69 { AbsMult(x, y); }
70
71 /// @deprecated Use AbsMultTranspose() instead.
72 MFEM_DEPRECATED void MultTransposeUnsigned(const Vector &x, Vector &y) const
73 { AbsMultTranspose(x, y); }
74
75 /// Compute MultTranspose by setting (rather than adding) element
76 /// contributions; this is a left inverse of the Mult() operation
77 void MultLeftInverse(const Vector &x, Vector &y) const;
78
79 /// @brief Fills the E-vector y with `boolean` values 0.0 and 1.0 such that each
80 /// each entry of the L-vector is uniquely represented in `y`.
81 /** This means, the sum of the E-vector `y` is equal to the sum of the
82 corresponding L-vector filled with ones. The boolean mask is required to
83 emulate SetSubVector and its transpose on GPUs. This method is running on
84 the host, since the `processed` array requires a large shared memory. */
85 void BooleanMask(Vector& y) const;
86
87 /// Fill a Sparse Matrix with Element Matrices.
88 void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const;
89
90 /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
91 given by this ElementRestriction. */
92 int FillI(SparseMatrix &mat) const;
93 /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
94 pattern given by this ElementRestriction, and the values of ea_data. */
95 void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const;
96 /// @private Not part of the public interface (device kernel limitation).
97 ///
98 /// Performs either MultTranspose or AddMultTranspose depending on the
99 /// boolean template parameter @a ADD.
100 template <bool ADD> void TAddMultTranspose(const Vector &x, Vector &y) const;
101
102 /// @name Low-level access to the underlying element-dof mappings
103 ///@{
104 const Array<int> &GatherMap() const { return gather_map; }
105 const Array<int> &Indices() const { return indices; }
106 const Array<int> &Offsets() const { return offsets; }
107 ///@}
108};
109
110/// Operator that converts L2 FiniteElementSpace L-vectors to E-vectors.
111/** Objects of this type are typically created and owned by FiniteElementSpace
112 objects, see FiniteElementSpace::GetElementRestriction(). L-vectors
113 corresponding to grid functions in L2 finite element spaces differ from
114 E-vectors only in the ordering of the degrees of freedom. */
116{
117 const int ne;
118 const int vdim;
119 const bool byvdim;
120 const int ndof;
121 const int ndofs;
122public:
124 void Mult(const Vector &x, Vector &y) const override;
125 void MultTranspose(const Vector &x, Vector &y) const override;
126 void AddMultTranspose(const Vector &x, Vector &y,
127 const real_t a = 1.0) const override;
128 /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
129 given by this ElementRestriction. */
130 void FillI(SparseMatrix &mat) const;
131 /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
132 pattern given by this L2FaceRestriction, and the values of ea_data. */
133 void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const;
134 /// @private Not part of the public interface (device kernel limitation).
135 ///
136 /// Performs either MultTranspose or AddMultTranspose depending on the
137 /// boolean template parameter @a ADD.
138 template <bool ADD> void TAddMultTranspose(const Vector &x, Vector &y) const;
139};
140
141/** An enum type to specify if only e1 value is requested (SingleValued) or both
142 e1 and e2 (DoubleValued). */
144
145/** @brief Base class for operators that extracts Face degrees of freedom.
146
147 In order to compute quantities on the faces of a mesh, it is often useful to
148 extract the degrees of freedom on the faces of the elements. This class
149 provides an interface for such operations.
150
151 If the FiniteElementSpace is ordered by Ordering::byVDIM, then the expected
152 format for the L-vector is (vdim x ndofs), otherwise if Ordering::byNODES
153 the expected format is (ndofs x vdim), where ndofs is the total number of
154 degrees of freedom.
155 Since FiniteElementSpace can either be continuous or discontinuous, the
156 degrees of freedom on a face can either be single valued or double valued,
157 this is what we refer to as the multiplicity and is represented by the
158 L2FaceValues enum type.
159 The format of the output face E-vector of degrees of freedom is
160 (face_dofs x vdim x multiplicity x nfaces), where face_dofs is the number of
161 degrees of freedom on each face, and nfaces the number of faces of the
162 requested FaceType (see FiniteElementSpace::GetNFbyType).
163
164 @note Objects of this type are typically created and owned by
165 FiniteElementSpace objects, see FiniteElementSpace::GetFaceRestriction(). */
167{
168public:
170
171 FaceRestriction(int h, int w): Operator(h, w) { }
172
173 virtual ~FaceRestriction() { }
174
175 /** @brief Extract the face degrees of freedom from @a x into @a y.
176
177 @param[in] x The L-vector of degrees of freedom.
178 @param[out] y The degrees of freedom on the face, corresponding to a face
179 E-vector.
180 */
181 void Mult(const Vector &x, Vector &y) const override = 0;
182
183 /** @brief Add the face degrees of freedom @a x to the element degrees of
184 freedom @a y.
185
186 @param[in] x The face degrees of freedom on the face.
187 @param[in,out] y The L-vector of degrees of freedom to which we add the
188 face degrees of freedom.
189 @param[in] a Scalar coefficient for addition.
190 */
191 void AddMultTranspose(const Vector &x, Vector &y,
192 const real_t a = 1.0) const override = 0;
193
194 /** @brief Add the face degrees of freedom @a x to the element degrees of
195 freedom @a y ignoring the signs from DOF orientation. */
196 virtual void AddAbsMultTranspose(const Vector &x, Vector &y,
197 const real_t a = 1.0) const
198 {
199 AddMultTranspose(x, y, a);
200 }
201
202 /// @deprecated Use AddAbsMultTranspose() instead.
203 MFEM_DEPRECATED void AddMultTransposeUnsigned(const Vector &x, Vector &y,
204 const real_t a = 1.0) const
205 {
206 AddAbsMultTranspose(x, y, a);
207 }
208
209 /** @brief Add the face degrees of freedom @a x to the element degrees of
210 freedom @a y. Perform the same computation as AddMultTranspose, but
211 @a x is invalid after calling this method.
212
213 @param[in,out] x The face degrees of freedom on the face.
214 @param[in,out] y The L-vector of degrees of freedom to which we add the
215 face degrees of freedom.
216
217 @note This method is an optimization of AddMultTranspose where the @a x
218 Vector is used and modified to avoid memory allocation and memcpy.
219 */
220 virtual void AddMultTransposeInPlace(Vector &x, Vector &y) const
221 {
222 AddMultTranspose(x, y);
223 }
224
225 /** @brief Set the face degrees of freedom in the element degrees of freedom
226 @a y to the values given in @a x.
227
228 @param[in] x The face degrees of freedom on the face.
229 @param[in,out] y The L-vector of degrees of freedom to which we add the
230 face degrees of freedom.
231 */
232 void MultTranspose(const Vector &x, Vector &y) const override
233 {
234 y = 0.0;
235 AddMultTranspose(x, y);
236 }
237
238 void AbsMultTranspose(const Vector &x, Vector &y) const override
239 {
240 y = 0.0;
242 }
243
244 /** @brief For each face, sets @a y to the partial derivative of @a x with
245 respect to the reference coordinate whose direction is
246 perpendicular to the face on the reference element.
247
248 @details This is not the normal derivative in physical coordinates, but can
249 be mapped to the physical normal derivative using the element
250 Jacobian and the tangential derivatives (in reference coordinates)
251 which can be computed from the face values (provided by Mult).
252
253 Note that due to the polynomial degree of the element mapping, the
254 physical normal derivative may be a higher degree polynomial than
255 the restriction of the values to the face. However, the normal
256 derivative in reference coordinates has degree-1, and therefore can
257 be exactly represented with the degrees of freedom of a face
258 E-vector.
259
260 @param[in] x The L-vector degrees of freedom.
261 @param[in,out] y The reference normal derivative degrees of freedom. Is
262 E-vector like.
263 */
264 virtual void NormalDerivativeMult(const Vector &x, Vector &y) const
265 {
266 MFEM_ABORT("Not implemented for this restriction operator.");
267 }
268
269 /** @brief Add the face reference-normal derivative degrees of freedom in @a
270 x to the element degrees of freedom in @a y.
271
272 @details see NormalDerivativeMult.
273
274 @param[in] x The degrees of freedom of the face reference-normal
275 derivative. Is E-vector like.
276 @param[in,out] y The L-vector degrees of freedom.
277 */
278 virtual void NormalDerivativeAddMultTranspose(const Vector &x, Vector &y) const
279 {
280 MFEM_ABORT("Not implemented for this restriction operator.");
281 }
282
283 /// @brief Low-level access to the underlying gather map.
284 virtual const Array<int> &GatherMap() const
285 {
286 MFEM_ABORT("Not implemented for this restriction operator.");
287 }
288};
289
290/// @brief Operator that extracts face degrees of freedom for H1, ND, or RT
291/// FiniteElementSpaces.
292///
293/// Objects of this type are typically created and owned by FiniteElementSpace
294/// objects, see FiniteElementSpace::GetFaceRestriction().
296{
297protected:
299 const int nf; // Number of faces of the requested type
300 const int vdim;
301 const bool byvdim;
302 const int face_dofs; // Number of dofs on each face
303 const int elem_dofs; // Number of dofs in each element
304 const int nfdofs; // Total number of face E-vector dofs
305 const int ndofs; // Total number of dofs
306 Array<int> scatter_indices; // Scattering indices for element 1 on each face
307 Array<int> gather_offsets; // offsets for the gathering indices of each dof
308 Array<int> gather_indices; // gathering indices for each dof
309 Array<int> vol_dof_map; // mapping from lexicographic to native ordering
310
311 /** @brief Construct a ConformingFaceRestriction.
312
313 @param[in] fes The FiniteElementSpace on which this operates
314 @param[in] f_ordering Request a specific face dof ordering
315 @param[in] type Request internal or boundary faces dofs
316 @param[in] build Request the NCL2FaceRestriction to compute the
317 scatter/gather indices. False should only be used
318 when inheriting from ConformingFaceRestriction.
319 */
321 const ElementDofOrdering f_ordering,
322 const FaceType type,
323 bool build);
324public:
325 /** @brief Construct a ConformingFaceRestriction.
326
327 @param[in] fes The FiniteElementSpace on which this operates
328 @param[in] f_ordering Request a specific face dof ordering
329 @param[in] type Request internal or boundary faces dofs */
331 const ElementDofOrdering f_ordering,
332 const FaceType type);
333
334 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
335 face E-Vector.
336
337 @param[in] x The L-vector degrees of freedom.
338 @param[out] y The face E-Vector degrees of freedom with the given format:
339 face_dofs x vdim x nf
340 where nf is the number of interior or boundary faces
341 requested by @a type in the constructor.
342 The face_dofs are ordered according to the given
343 ElementDofOrdering. */
344 void Mult(const Vector &x, Vector &y) const override
345 { MultInternal(x, y); }
346
347 /// Compute Mult without applying signs based on DOF orientations.
348 void AbsMult(const Vector &x, Vector &y) const override
349 { MultInternal(x, y, true); }
350
351 /// @deprecated Use AbsMult() instead.
352 MFEM_DEPRECATED void MultUnsigned(const Vector &x, Vector &y) const
353 { AbsMult(x, y); }
354
356
357 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
358 L-Vector.
359
360 @param[in] x The face E-Vector degrees of freedom with the given format:
361 face_dofs x vdim x nf
362 where nf is the number of interior or boundary faces
363 requested by @a type in the constructor.
364 The face_dofs should be ordered according to the given
365 ElementDofOrdering
366 @param[in,out] y The L-vector degrees of freedom.
367 @param[in] a Scalar coefficient for addition. */
368 void AddMultTranspose(const Vector &x, Vector &y,
369 const real_t a = 1.0) const override;
370
371 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
372 L-Vector @b not taking into account signs from DOF orientations.
373
374 @sa AddMultTranspose(). */
375 void AddAbsMultTranspose(const Vector &x, Vector &y,
376 const real_t a = 1.0) const override;
377
378 /// @deprecated Use AddAbsMultTranspose() instead.
379 MFEM_DEPRECATED void AddMultTransposeUnsigned(const Vector &x, Vector &y) const
380 {
382 }
383
384 void AbsMultTranspose(const Vector &x, Vector &y) const override
385 {
386 y = 0.0;
388 }
389
390private:
391 /** @brief Compute the scatter indices: L-vector to E-vector, and the offsets
392 for the gathering: E-vector to L-vector.
393
394 @param[in] f_ordering Request a specific face dof ordering.
395 @param[in] type Request internal or boundary faces dofs.
396 */
397 void ComputeScatterIndicesAndOffsets(const ElementDofOrdering f_ordering,
398 const FaceType type);
399
400 /** @brief Compute the gather indices: E-vector to L-vector.
401
402 Note: Requires the gather offsets to be computed.
403
404 @param[in] f_ordering Request a specific face dof ordering.
405 @param[in] type Request internal or boundary faces dofs.
406 */
407 void ComputeGatherIndices(const ElementDofOrdering f_ordering,
408 const FaceType type);
409
410protected:
411 mutable Array<int> face_map; // Used in the computation of GetFaceDofs
412
413 /** @brief Verify that ConformingFaceRestriction is built from a supported
414 finite element space.
415
416 @param[in] f_ordering The requested face dof ordering.
417 */
418 void CheckFESpace(const ElementDofOrdering f_ordering);
419
420 /** @brief Set the scattering indices of elem1, and increment the offsets for
421 the face described by the @a face.
422
423 @param[in] face The face information of the current face.
424 @param[in] face_index The interior/boundary face index.
425 @param[in] f_ordering Request a specific face dof ordering.
426 */
428 const int face_index,
429 const ElementDofOrdering f_ordering);
430
431 /** @brief Set the gathering indices of elem1 for the interior face described
432 by the @a face.
433
434 @param[in] face The face information of the current face.
435 @param[in] face_index The interior/boundary face index.
436 @param[in] f_ordering Request a specific face dof ordering.
437 */
439 const int face_index,
440 const ElementDofOrdering f_ordering);
441
442public:
443 // This method needs to be public due to 'nvcc' restriction.
444 void MultInternal(const Vector &x, Vector &y,
445 const bool useAbs = false) const;
446};
447
448/// @brief Alias for ConformingFaceRestriction, for backwards compatibility and
449/// as base class for ParNCH1FaceRestriction.
451
452/// Operator that extracts Face degrees of freedom for L2 spaces.
453/** Objects of this type are typically created and owned by FiniteElementSpace
454 objects, see FiniteElementSpace::GetFaceRestriction(). */
456{
457protected:
460 const int nf; // Number of faces of the requested type
461 const int ne; // Number of elements
462 const int vdim; // vdim
463 const bool byvdim;
464 const int face_dofs; // Number of dofs on each face
465 const int elem_dofs; // Number of dofs in each element
466 const int nfdofs; // Total number of dofs on the faces
467 const int ndofs; // Total number of dofs
470 Array<int> scatter_indices1; // Scattering indices for element 1 on each face
471 Array<int> scatter_indices2; // Scattering indices for element 2 on each face
472 Array<int> gather_offsets; // offsets for the gathering indices of each dof
473 Array<int> gather_indices; // gathering indices for each dof
474 mutable std::unique_ptr<L2NormalDerivativeFaceRestriction> normal_deriv_restr;
475
476 /** @brief Constructs an L2FaceRestriction.
477
478 @param[in] fes The FiniteElementSpace on which this operates
479 @param[in] f_ordering Request a specific face dof ordering
480 @param[in] type Request internal or boundary faces dofs
481 @param[in] m Request the face dofs for elem1, or both elem1 and
482 elem2
483 @param[in] build Request the NCL2FaceRestriction to compute the
484 scatter/gather indices. False should only be used
485 when inheriting from L2FaceRestriction.
486 */
488 const ElementDofOrdering f_ordering,
489 const FaceType type,
490 const L2FaceValues m,
491 bool build);
492
493public:
494 /** @brief Constructs an L2FaceRestriction.
495
496 @param[in] fes The FiniteElementSpace on which this operates
497 @param[in] f_ordering Request a specific face dof ordering
498 @param[in] type Request internal or boundary faces dofs
499 @param[in] m Request the face dofs for elem1, or both elem1 and
500 elem2 */
502 const ElementDofOrdering f_ordering,
503 const FaceType type,
505
506 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
507 face E-Vector.
508
509 @param[in] x The L-vector degrees of freedom.
510 @param[out] y The face E-Vector degrees of freedom with the given format:
511 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf)
512 if L2FacesValues::SingleValued (face_dofs x vdim x nf)
513 where nf is the number of interior or boundary faces
514 requested by @a type in the constructor.
515 The face_dofs are ordered according to the given
516 ElementDofOrdering. */
517 void Mult(const Vector &x, Vector &y) const override;
518
520
521 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
522 L-Vector.
523
524 @param[in] x The face E-Vector degrees of freedom with the given format:
525 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf)
526 if L2FacesValues::SingleValued (face_dofs x vdim x nf)
527 where nf is the number of interior or boundary faces
528 requested by @a type in the constructor.
529 The face_dofs should be ordered according to the given
530 ElementDofOrdering
531 @param[in,out] y The L-vector degrees of freedom.
532 @param[in] a Scalar coefficient for addition. */
533 void AddMultTranspose(const Vector &x, Vector &y,
534 const real_t a = 1.0) const override;
535
536 /** @brief Fill the I array of SparseMatrix corresponding to the sparsity
537 pattern given by this L2FaceRestriction.
538
539 @param[in,out] mat The sparse matrix for which we want to initialize the
540 row offsets.
541 @param[in] keep_nbr_block When set to true the SparseMatrix will
542 include the rows (in addition to the columns)
543 corresponding to face-neighbor dofs. The
544 default behavior is to disregard those rows. */
545 virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block = false) const;
546
547 /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
548 the sparsity pattern given by this L2FaceRestriction, and the values of
549 fea_data.
550
551 @param[in] fea_data The dense matrices representing the local operators
552 on each face. The format is:
553 face_dofs x face_dofs x 2 x nf
554 On each face the first local matrix corresponds to
555 the contribution of elem1 on elem2, and the second to
556 the contribution of elem2 on elem1.
557 @param[in,out] mat The sparse matrix that is getting filled.
558 @param[in] keep_nbr_block When set to true the SparseMatrix will
559 include the rows (in addition to the columns)
560 corresponding to face-neighbor dofs. The
561 default behavior is to disregard those rows. */
562 virtual void FillJAndData(const Vector &fea_data,
563 SparseMatrix &mat,
564 const bool keep_nbr_block = false) const;
565
566 /** @brief This methods adds the DG face matrices to the element matrices.
567
568 @param[in] fea_data The dense matrices representing the local operators
569 on each face. The format is:
570 face_dofs x face_dofs x 2 x nf
571 On each face the first and second local matrices
572 correspond to the contributions of elem1 and elem2 on
573 themselves respectively.
574 @param[in,out] ea_data The dense matrices representing the element local
575 contributions for each element to which will be
576 added the face contributions.
577 The format is: dofs x dofs x ne, where dofs is the
578 number of dofs per element and ne the number of
579 elements. */
580 virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data,
581 Vector &ea_data) const;
582
583 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
584 face E-Vector.
585
586 @param[in] x The L-vector degrees of freedom.
587 @param[out] y The face E-Vector degrees of freedom with the given format:
588 (face_dofs x vdim x 2 x nf) where nf is the number of
589 interior or boundary faces requested by @a type in the
590 constructor. The face_dofs are ordered according to the
591 given ElementDofOrdering. */
592 void NormalDerivativeMult(const Vector &x, Vector &y) const override;
593
594 /** @brief Add the face reference-normal derivative degrees of freedom in @a
595 x to the element degrees of freedom in @a y.
596
597 @details see NormalDerivativeMult.
598
599 @param[in] x The degrees of freedom of the face reference-normal
600 derivative. Is E-vector like.
601 @param[in,out] y The L-vector degrees of freedom.
602 */
604 Vector &y) const override;
605private:
606 /** @brief Compute the scatter indices: L-vector to E-vector, and the offsets
607 for the gathering: E-vector to L-vector.
608 */
609 void ComputeScatterIndicesAndOffsets();
610
611 /** @brief Compute the gather indices: E-vector to L-vector.
612
613 Note: Requires the gather offsets to be computed.
614 */
615 void ComputeGatherIndices();
616
617 /// Create the internal normal derivative restriction operator if needed.
618 void EnsureNormalDerivativeRestriction() const;
619
620protected:
621 mutable Array<int> face_map; // Used in the computation of GetFaceDofs
622
623 /** @brief Verify that L2FaceRestriction is built from an L2 FESpace.
624 */
625 void CheckFESpace();
626
627 /** @brief Set the scattering indices of elem1, and increment the offsets for
628 the face described by the @a face. The ordering of the face dofs of elem1
629 is lexicographic relative to elem1.
630
631 @param[in] face The face information of the current face.
632 @param[in] face_index The interior/boundary face index.
633 */
635 const int face_index);
636
637 /** @brief Permute and set the scattering indices of elem2, and increment the
638 offsets for the face described by the @a face. The permutation orders the
639 dofs of elem2 lexicographically as the ones of elem1.
640
641 @param[in] face The face information of the current face.
642 @param[in] face_index The interior/boundary face index.
643 */
645 const int face_index);
646
647 /** @brief Permute and set the scattering indices of elem2 for the shared
648 face described by the @a face. The permutation orders the dofs of elem2 as
649 the ones of elem1.
650
651 @param[in] face The face information of the current face.
652 @param[in] face_index The interior/boundary face index.
653 */
655 const Mesh::FaceInformation &face,
656 const int face_index);
657
658 /** @brief Set the scattering indices of elem2 for the boundary face
659 described by the @a face.
660
661 @param[in] face The face information of the current face.
662 @param[in] face_index The interior/boundary face index.
663 */
665 const int face_index);
666
667 /** @brief Set the gathering indices of elem1 for the interior face described
668 by the @a face.
669
670 Note: This function modifies the offsets.
671
672 @param[in] face The face information of the current face.
673 @param[in] face_index The interior/boundary face index.
674 */
676 const int face_index);
677
678 /** @brief Permute and set the gathering indices of elem2 for the interior
679 face described by the @a face. The permutation orders the dofs of elem2 as
680 the ones of elem1.
681
682 Note: This function modifies the offsets.
683
684 @param[in] face The face information of the current face.
685 @param[in] face_index The interior/boundary face index.
686 */
688 const int face_index);
689
690public:
691 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
692 face E-Vector. Should only be used with conforming faces and when:
693 m == L2FacesValues::SingleValued
694
695 @param[in] x The L-vector degrees of freedom.
696 @param[out] y The face E-Vector degrees of freedom with the given format:
697 face_dofs x vdim x nf
698 where nf is the number of interior or boundary faces
699 requested by @a type in the constructor.
700 The face_dofs are ordered according to the given
701 ElementDofOrdering. */
702 void SingleValuedConformingMult(const Vector& x, Vector& y) const;
703
704 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
705 face E-Vector. Should only be used with conforming faces and when:
706 m == L2FacesValues::DoubleValued
707
708 @param[in] x The L-vector degrees of freedom.
709 @param[out] y The face E-Vector degrees of freedom with the given format:
710 face_dofs x vdim x 2 x nf
711 where nf is the number of interior or boundary faces
712 requested by @a type in the constructor.
713 The face_dofs are ordered according to the given
714 ElementDofOrdering. */
715 virtual void DoubleValuedConformingMult(const Vector& x, Vector& y) const;
716
717 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
718 L-Vector. Should only be used with conforming faces and when:
719 m == L2FacesValues::SingleValued
720
721 @param[in] x The face E-Vector degrees of freedom with the given format:
722 face_dofs x vdim x nf
723 where nf is the number of interior or boundary faces
724 requested by @a type in the constructor.
725 The face_dofs should be ordered according to the given
726 ElementDofOrdering
727 @param[in,out] y The L-vector degrees of freedom. */
729
730 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
731 L-Vector. Should only be used with conforming faces and when:
732 m == L2FacesValues::DoubleValued
733
734 @param[in] x The face E-Vector degrees of freedom with the given format:
735 face_dofs x vdim x 2 x nf
736 where nf is the number of interior or boundary faces
737 requested by @a type in the constructor.
738 The face_dofs should be ordered according to the given
739 ElementDofOrdering
740 @param[in,out] y The L-vector degrees of freedom. */
742};
743
744/** This struct stores which side is the master nonconforming side and the
745 index of the interpolator, see InterpolationManager class below. */
747{
748 uint32_t is_non_conforming : 1;
749 uint32_t master_side : 1;
750 uint32_t index : 30;
751
752 // default constructor, create a conforming face with index 0.
753 InterpConfig() = default;
754
755 // Non-conforming face
756 InterpConfig(int master_side, int nc_index)
758 { }
759
760 InterpConfig(const InterpConfig&) = default;
761
762 InterpConfig &operator=(const InterpConfig &rhs) = default;
763};
764
765/** This struct stores which side is the master nonconforming side and the
766 index of the interpolator, see InterpolationManager class below. */
768{
770 uint32_t is_non_conforming : 1;
771 uint32_t master_side : 1;
772 uint32_t index : 30;
773
774 // default constructor.
775 NCInterpConfig() = default;
776
777 // Non-conforming face
778 NCInterpConfig(int face_index, int master_side, int nc_index)
782 index(nc_index)
783 { }
784
785 // Non-conforming face
792
794
796};
797
798/** @brief This class manages the storage and computation of the interpolations
799 from master (coarse) face to slave (fine) face.
800*/
802{
803protected:
806 Array<InterpConfig> interp_config; // interpolator index for each face
807 Array<NCInterpConfig> nc_interp_config; // interpolator index for each ncface
808 Vector interpolators; // face_dofs x face_dofs x num_interpolators
809 int nc_cpt; // Counter for interpolators, and used as index.
810
811 /** The interpolators are associated to a key of containing the address of
812 PointMatrix and a local face identifier. */
813 using Key = std::pair<const DenseMatrix*,int>;
814 /// The temporary map used to store the different interpolators.
815 using Map = std::map<Key, std::pair<int,const DenseMatrix*>>;
816 Map interp_map; // The temporary map that stores the interpolators.
817
818public:
820
821 /** @brief main constructor.
822
823 @param[in] fes The FiniteElementSpace on which this operates
824 @param[in] ordering Request a specific element ordering.
825 @param[in] type Request internal or boundary faces dofs
826 */
829 FaceType type);
830
831 /** @brief Register the face with @a face and index @a face_index as a
832 conforming face for the interpolation of the degrees of freedom.
833
834 @param[in] face The face information of the current face.
835 @param[in] face_index The interior/boundary face index.
836 */
838 int face_index);
839
840 /** @brief Register the face with @a face and index @a face_index as a
841 conforming face for the interpolation of the degrees of freedom.
842
843 @param[in] face The face information of the current face.
844 @param[in] face_index The interior/boundary face index.
845 */
847 int face_index);
848
849 /** @brief Transform the interpolation matrix map into a contiguous memory
850 structure. */
852
854
855 /// @brief Return the total number of interpolators.
857 {
858 return nc_cpt;
859 }
860
861 /** @brief Return an mfem::Vector containing the interpolators in the
862 following format: face_dofs x face_dofs x num_interpolators. */
864 {
865 return interpolators;
866 }
867
868 /** @brief Return an array containing the interpolation configuration for
869 each face registered with RegisterFaceConformingInterpolation and
870 RegisterFaceCoarseToFineInterpolation. */
872 {
873 return interp_config;
874 }
875
876 /** @brief Return an array containing the interpolation configuration for
877 each face registered with RegisterFaceConformingInterpolation and
878 RegisterFaceCoarseToFineInterpolation. */
880 {
881 return nc_interp_config;
882 }
883
884private:
885 /** @brief Returns the interpolation operator from a master (coarse) face to
886 a slave (fine) face.
887
888 @param[in] face The face information of the current face.
889 @param[in] ptMat The PointMatrix describing the position and orientation
890 of the fine face in the coarse face. This PointMatrix is
891 usually obtained from the mesh through the method
892 GetNCFacesPtMat.
893 @param[in] ordering Request a specific element ordering.
894 @return The dense matrix corresponding to the interpolation of the face
895 degrees of freedom of the master (coarse) face to the slave
896 (fine) face. */
897 const DenseMatrix* GetCoarseToFineInterpolation(
898 const Mesh::FaceInformation &face,
899 const DenseMatrix* ptMat);
900};
901
902/** @brief Operator that extracts face degrees of freedom for L2 nonconforming
903 spaces.
904
905 In order to support face restrictions on nonconforming meshes, this
906 operator interpolates master (coarse) face degrees of freedom onto the
907 slave (fine) face. This allows face integrators to treat nonconforming
908 faces just as regular conforming faces. */
910{
911protected:
914
915 /** @brief Constructs an NCL2FaceRestriction, this is a specialization of a
916 L2FaceRestriction for nonconforming meshes.
917
918 @param[in] fes The FiniteElementSpace on which this operates
919 @param[in] f_ordering Request a specific face dof ordering
920 @param[in] type Request internal or boundary faces dofs
921 @param[in] m Request the face dofs for elem1, or both elem1 and
922 elem2
923 @param[in] build Request the NCL2FaceRestriction to compute the
924 scatter/gather indices. False should only be used
925 when inheriting from NCL2FaceRestriction.
926 */
928 const ElementDofOrdering f_ordering,
929 const FaceType type,
930 const L2FaceValues m,
931 bool build);
932public:
933 /** @brief Constructs an NCL2FaceRestriction, this is a specialization of a
934 L2FaceRestriction for nonconforming meshes.
935
936 @param[in] fes The FiniteElementSpace on which this operates
937 @param[in] f_ordering Request a specific face dof ordering
938 @param[in] type Request internal or boundary faces dofs
939 @param[in] m Request the face dofs for elem1, or both elem1 and
940 elem2
941 */
943 const ElementDofOrdering f_ordering,
944 const FaceType type,
946
947 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
948 face E-Vector.
949
950 @param[in] x The L-vector degrees of freedom.
951 @param[out] y The face E-Vector degrees of freedom with the given format:
952 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
953 if L2FacesValues::SingleValued (face_dofs x vdim x nf),
954 where nf is the number of interior or boundary faces
955 requested by @a type in the constructor.
956 The face_dofs are ordered according to the given
957 ElementDofOrdering. */
958 void Mult(const Vector &x, Vector &y) const override;
959
960 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
961 L-Vector.
962
963 @param[in] x The face E-Vector degrees of freedom with the given format:
964 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
965 if L2FacesValues::SingleValued (face_dofs x vdim x nf),
966 where nf is the number of interior or boundary faces
967 requested by @a type in the constructor.
968 The face_dofs should be ordered according to the given
969 ElementDofOrdering
970 @param[in,out] y The L-vector degrees of freedom.
971 @param[in] a Scalar coefficient for addition. */
972 void AddMultTranspose(const Vector &x, Vector &y,
973 const real_t a = 1.0) const override;
974
975 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
976 L-Vector.
977
978 @param[in,out] x The face E-Vector degrees of freedom with the given format:
979 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
980 if L2FacesValues::SingleValued (face_dofs x vdim x nf),
981 where nf is the number of interior or boundary faces
982 requested by @a type in the constructor.
983 The face_dofs should be ordered according to the given
984 ElementDofOrdering
985 @param[in,out] y The L-vector degrees of freedom.
986
987 @note This method is an optimization of AddMultTranspose where the @a x
988 Vector is used and modified to avoid memory allocation and memcpy. */
989 void AddMultTransposeInPlace(Vector &x, Vector &y) const override;
990
991 /** @brief Fill the I array of SparseMatrix corresponding to the sparsity
992 pattern given by this NCL2FaceRestriction.
993
994 @param[in,out] mat The sparse matrix for which we want to initialize the
995 row offsets.
996 @param[in] keep_nbr_block When set to true the SparseMatrix will
997 include the rows (in addition to the columns)
998 corresponding to face-neighbor dofs. The
999 default behavior is to disregard those rows.
1000
1001 @warning This method is not implemented yet. */
1002 void FillI(SparseMatrix &mat,
1003 const bool keep_nbr_block = false) const override;
1004
1005 /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
1006 the sparsity pattern given by this NCL2FaceRestriction, and the values of
1007 ea_data.
1008
1009 @param[in] fea_data The dense matrices representing the local operators
1010 on each face. The format is:
1011 face_dofs x face_dofs x 2 x nf.
1012 On each face the first local matrix corresponds to
1013 the contribution of elem1 on elem2, and the second to
1014 the contribution of elem2 on elem1.
1015 @param[in,out] mat The sparse matrix that is getting filled.
1016 @param[in] keep_nbr_block When set to true the SparseMatrix will
1017 include the rows (in addition to the columns)
1018 corresponding to face-neighbor dofs. The
1019 default behavior is to disregard those rows.
1020
1021 @warning This method is not implemented yet. */
1022 void FillJAndData(const Vector &fea_data,
1023 SparseMatrix &mat,
1024 const bool keep_nbr_block = false) const override;
1025
1026 /** @brief This methods adds the DG face matrices to the element matrices.
1027
1028 @param[in] fea_data The dense matrices representing the local operators
1029 on each face. The format is:
1030 face_dofs x face_dofs x 2 x nf.
1031 On each face the first and second local matrices
1032 correspond to the contributions of elem1 and elem2 on
1033 themselves respectively.
1034 @param[in,out] ea_data The dense matrices representing the element local
1035 contributions for each element to which will be
1036 added the face contributions.
1037 The format is: dofs x dofs x ne, where dofs is the
1038 number of dofs per element and ne the number of
1039 elements.
1040
1041 @warning This method is not implemented yet. */
1042 void AddFaceMatricesToElementMatrices(const Vector &fea_data,
1043 Vector &ea_data) const override;
1044
1045private:
1046 /** @brief Compute the scatter indices: L-vector to E-vector, the offsets
1047 for the gathering: E-vector to L-vector, and the interpolators from
1048 coarse to fine face for master non-comforming faces.
1049 */
1050 void ComputeScatterIndicesAndOffsets();
1051
1052 /** @brief Compute the gather indices: E-vector to L-vector.
1053
1054 Note: Requires the gather offsets to be computed.
1055 */
1056 void ComputeGatherIndices();
1057
1058public:
1059 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
1060 face E-Vector. Should only be used with nonconforming faces and when:
1061 L2FaceValues m == L2FaceValues::DoubleValued
1062
1063 @param[in] x The L-vector degrees of freedom.
1064 @param[out] y The face E-Vector degrees of freedom with the given format:
1065 (face_dofs x vdim x 2 x nf),
1066 where nf is the number of interior or boundary faces
1067 requested by @a type in the constructor.
1068 The face_dofs are ordered according to the given
1069 ElementDofOrdering. */
1070 virtual void DoubleValuedNonconformingMult(const Vector& x, Vector& y) const;
1071
1072 /** @brief Apply a change of basis from coarse element basis to fine element
1073 basis for the coarse face dofs.
1074
1075 @param[in,out] x The dofs vector that needs coarse dofs to be express in
1076 term of the fine basis.
1077 */
1079
1080 /** @brief Apply a change of basis from fine element basis to coarse element
1081 basis for the coarse face dofs. Should only be used when:
1082 L2FaceValues m == L2FaceValues::SingleValued
1083
1084 @param[in] x The dofs vector that needs coarse dofs to be express in term
1085 of the coarse basis, the result is stored in x_interp.
1086 */
1088
1089 /** @brief Apply a change of basis from fine element basis to coarse element
1090 basis for the coarse face dofs. Should only be used when:
1091 L2FaceValues m == L2FaceValues::SingleValued
1092
1093 @param[in,out] x The dofs vector that needs coarse dofs to be express in
1094 term of the coarse basis, the result is stored in x.
1095 */
1097
1098 /** @brief Apply a change of basis from fine element basis to coarse element
1099 basis for the coarse face dofs. Should only be used when:
1100 L2FaceValues m == L2FaceValues::DoubleValued
1101
1102 @param[in] x The dofs vector that needs coarse dofs to be express in term
1103 of the coarse basis, the result is stored in x_interp.
1104 */
1106
1107 /** @brief Apply a change of basis from fine element basis to coarse element
1108 basis for the coarse face dofs. Should only be used when:
1109 L2FaceValues m == L2FaceValues::DoubleValued
1110
1111 @param[in,out] x The dofs vector that needs coarse dofs to be express in
1112 term of the coarse basis, the result is stored in
1113 x.
1114 */
1116};
1117
1118/// Operator that extracts face degrees of freedom for L2 interface spaces.
1119/** Objects of this type are typically created and owned by FiniteElementSpace
1120 objects, see FiniteElementSpace::GetFaceRestriction(). */
1122{
1123protected:
1124 const FiniteElementSpace &fes; ///< The finite element space
1125 const ElementDofOrdering ordering; ///< Requested ordering
1126 const FaceType type; ///< Face type (interior or boundary)
1127 const int nfaces; ///< Number of faces of the requested type
1128 const int vdim; ///< vdim of the space
1129 const bool byvdim; ///< DOF ordering (by nodes or by vdim)
1130 const int face_dofs; ///< Number of dofs on each face
1131 const int nfdofs; ///< Total number of dofs on the faces (E-vector size)
1132 const int ndofs; ///< Number of dofs in the space (L-vector size)
1133 Array<int> gather_map; ///< Gather map
1134
1135public:
1136 /** @brief Constructs an L2InterfaceFaceRestriction.
1137
1138 @param[in] fes_ The FiniteElementSpace on which this operates
1139 @param[in] ordering_ Request a specific face dof ordering
1140 @param[in] type_ Request internal or boundary faces dofs */
1142 const ElementDofOrdering ordering_,
1143 const FaceType type_);
1144
1145 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
1146 face E-Vector.
1147
1148 @param[in] x The L-vector degrees of freedom.
1149 @param[out] y The face E-Vector degrees of freedom with size (face_dofs,
1150 vdim, nf), where nf is the number of interior or boundary
1151 faces requested by @a type in the constructor. The
1152 face_dofs are ordered according to the given
1153 ElementDofOrdering. */
1154 void Mult(const Vector &x, Vector &y) const override;
1155
1157
1158 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
1159 L-Vector.
1160
1161 @param[in] x The face E-Vector degrees of freedom with size
1162 (face_dofs, vdim, nf), where nf is the number of
1163 interior or boundary faces requested by @a type in the
1164 constructor. The face_dofs should be ordered according
1165 to the given ElementDofOrdering
1166 @param[in,out] y The L-vector degrees of freedom.
1167 @param[in] a Scalar coefficient for addition. */
1168 void AddMultTranspose(const Vector &x, Vector &y,
1169 const real_t a = 1.0) const override;
1170
1171 const Array<int> &GatherMap() const override;
1172};
1173
1174/** @brief Convert a dof face index from Native ordering to lexicographic
1175 ordering for quads and hexes.
1176
1177 @param[in] dim The dimension of the element, 2 for quad, 3 for hex
1178 @param[in] face_id The local face identifier
1179 @param[in] size1d The 1D number of degrees of freedom for each dimension
1180 @param[in] index The native index on the face
1181 @return The lexicographic index on the face
1182*/
1183int ToLexOrdering(const int dim, const int face_id, const int size1d,
1184 const int index);
1185
1186/** @brief Compute the dof face index of elem2 corresponding to the given dof
1187 face index.
1188
1189 @param[in] dim The dimension of the element, 2 for quad, 3 for hex
1190 @param[in] face_id1 The local face identifier of elem1
1191 @param[in] face_id2 The local face identifier of elem2
1192 @param[in] orientation The orientation of elem2 relative to elem1 on the
1193 face
1194 @param[in] size1d The 1D number of degrees of freedom for each dimension
1195 @param[in] index The dof index on elem1
1196 @return The dof index on elem2 facing the dof on elem1
1197*/
1198int PermuteFaceL2(const int dim, const int face_id1,
1199 const int face_id2, const int orientation,
1200 const int size1d, const int index);
1201
1202/// @brief Return the face-neighbor data given the L-vector @a x.
1203///
1204/// If the input vector @a x is a ParGridFunction with non-empty face-neighbor
1205/// data, return an alias to ParGridFunction::FaceNbrData() (avoiding an
1206/// unneeded call to ParGridFunction::ExchangeFaceNbrData).
1207///
1208/// Otherwise, create a temporary ParGridFunction, exchange the face-neighbor
1209/// data, and return the resulting vector.
1210///
1211/// If @a fes is not a parallel space, or if @a ftype is not FaceType::Interior,
1212/// return an empty vector.
1214 const FiniteElementSpace &fes, const Vector &x, FaceType ftype);
1215
1216}
1217
1218#endif // MFEM_RESTRICTION
Operator that extracts face degrees of freedom for H1, ND, or RT FiniteElementSpaces.
ConformingFaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, bool build)
Construct a ConformingFaceRestriction.
void CheckFESpace(const ElementDofOrdering f_ordering)
Verify that ConformingFaceRestriction is built from a supported finite element space.
void SetFaceDofsGatherIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering f_ordering)
Set the gathering indices of elem1 for the interior face described by the face.
MFEM_DEPRECATED void AddMultTransposeUnsigned(const Vector &x, Vector &y) const
void MultInternal(const Vector &x, Vector &y, const bool useAbs=false) const
void AbsMult(const Vector &x, Vector &y) const override
Compute Mult without applying signs based on DOF orientations.
const FiniteElementSpace & fes
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void SetFaceDofsScatterIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering f_ordering)
Set the scattering indices of elem1, and increment the offsets for the face described by the face.
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 AbsMultTranspose(const Vector &x, Vector &y) const override
Action of the transpose absolute-value operator: y=|A|^t(x). The default behavior in class Operator i...
MFEM_DEPRECATED void MultUnsigned(const Vector &x, Vector &y) const
void AddAbsMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector not taking into account signs...
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Abstract base class that defines an interface for element restrictions.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override=0
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
Operator that converts FiniteElementSpace L-vectors to E-vectors.
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
const FiniteElementSpace & fes
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
const Array< int > & GatherMap() const
void AbsMult(const Vector &x, Vector &y) const override
Compute Mult without applying signs based on DOF orientations.
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
MFEM_DEPRECATED void MultTransposeUnsigned(const Vector &x, Vector &y) const
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
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 ...
const Array< int > & Offsets() const
void MultLeftInverse(const Vector &x, Vector &y) const
const Array< int > & Indices() const
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...
MFEM_DEPRECATED void MultUnsigned(const Vector &x, Vector &y) const
void AbsMultTranspose(const Vector &x, Vector &y) const override
Compute MultTranspose without applying signs based on DOF orientations.
int FillI(SparseMatrix &mat) const
Base class for operators that extracts Face degrees of freedom.
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 AbsMultTranspose(const Vector &x, Vector &y) const override
Action of the transpose absolute-value operator: y=|A|^t(x). The default behavior in class Operator i...
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override=0
Add the face degrees of freedom x to the element degrees of freedom y.
virtual void NormalDerivativeAddMultTranspose(const Vector &x, Vector &y) const
Add the face reference-normal derivative degrees of freedom in x to the element degrees of freedom in...
FaceRestriction(int h, int w)
virtual void NormalDerivativeMult(const Vector &x, Vector &y) const
For each face, sets y to the partial derivative of x with respect to the reference coordinate whose d...
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.
MFEM_DEPRECATED void AddMultTransposeUnsigned(const Vector &x, Vector &y, const real_t a=1.0) const
virtual void AddAbsMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const
Add the face degrees of freedom x to the element degrees of freedom y ignoring the signs from DOF ori...
virtual const Array< int > & GatherMap() const
Low-level access to the underlying gather map.
void Mult(const Vector &x, Vector &y) const override=0
Extract the face degrees of freedom from x into y.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
This class manages the storage and computation of the interpolations from master (coarse) face to sla...
std::map< Key, std::pair< int, const DenseMatrix * > > Map
The temporary map used to store the different interpolators.
std::pair< const DenseMatrix *, int > Key
const Array< NCInterpConfig > & GetNCFaceInterpConfig() const
Return an array containing the interpolation configuration for each face registered with RegisterFace...
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...
int GetNumInterpolators() const
Return the total number of interpolators.
Array< NCInterpConfig > nc_interp_config
const Vector & GetInterpolators() const
Return an mfem::Vector containing the interpolators in the following format: face_dofs x face_dofs x ...
void LinearizeInterpolatorMapIntoVector()
Transform the interpolation matrix map into a contiguous memory structure.
const FiniteElementSpace & fes
const ElementDofOrdering ordering
Array< InterpConfig > interp_config
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 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.
L2ElementRestriction(const FiniteElementSpace &)
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
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 ...
void FillI(SparseMatrix &mat) const
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 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....
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...
Array< int > scatter_indices2
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...
L2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an L2FaceRestriction.
void NormalDerivativeMult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void NormalDerivativeAddMultTranspose(const Vector &x, Vector &y) const override
Add the face reference-normal derivative degrees of freedom in x to the element degrees of freedom in...
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 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...
std::unique_ptr< L2NormalDerivativeFaceRestriction > normal_deriv_restr
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....
void CheckFESpace()
Verify that L2FaceRestriction is built from an L2 FESpace.
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...
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 ...
Array< int > scatter_indices1
const L2FaceValues m
const FiniteElementSpace & fes
void SetBoundaryDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem2 for the boundary face described by the face.
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 SetFaceDofsGatherIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the gathering indices of elem1 for the interior face described by the face.
const ElementDofOrdering ordering
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...
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...
virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
Operator that extracts face degrees of freedom for L2 interface spaces.
const int face_dofs
Number of dofs on each face.
const int ndofs
Number of dofs in the space (L-vector size)
const int nfdofs
Total number of dofs on the faces (E-vector size)
const FiniteElementSpace & fes
The finite element space.
L2InterfaceFaceRestriction(const FiniteElementSpace &fes_, const ElementDofOrdering ordering_, const FaceType type_)
Constructs an L2InterfaceFaceRestriction.
Array< int > gather_map
Gather map.
const int vdim
vdim of the space
const int nfaces
Number of faces of the requested type.
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
const Array< int > & GatherMap() const override
Low-level access to the underlying gather map.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
const ElementDofOrdering ordering
Requested ordering.
const FaceType type
Face type (interior or boundary)
const bool byvdim
DOF ordering (by nodes or by vdim)
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
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...
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 &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...
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
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....
InterpolationManager interpolations
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...
void DoubleValuedNonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs....
void DoubleValuedNonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs.
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 AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void DoubleValuedNonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs....
NCL2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an NCL2FaceRestriction, this is a specialization of a L2FaceRestriction for nonconforming ...
Abstract operator.
Definition operator.hpp:25
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:82
int dim
Definition ex24.cpp:53
int index(int i, int j, int nx, int ny)
Definition life.cpp:236
real_t a
Definition lissajous.cpp:41
Vector GetLVectorFaceNbrData(const FiniteElementSpace &fes, const Vector &x, FaceType ftype)
Return the face-neighbor data given the L-vector x.
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.
float real_t
Definition config.hpp:46
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.
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:47
FaceType
Definition mesh.hpp:49
InterpConfig & operator=(const InterpConfig &rhs)=default
InterpConfig(int master_side, int nc_index)
InterpConfig(const InterpConfig &)=default
InterpConfig()=default
This structure is used as a human readable output format that deciphers the information contained in ...
Definition mesh.hpp:2081
NCInterpConfig & operator=(const NCInterpConfig &rhs)=default
NCInterpConfig(const NCInterpConfig &)=default
NCInterpConfig(int face_index, InterpConfig &config)
NCInterpConfig()=default
NCInterpConfig(int face_index, int master_side, int nc_index)