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