MFEM  v4.6.0
Finite element discretization library
prestriction.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_PRESTRICTION
13 #define MFEM_PRESTRICTION
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include "restriction.hpp"
20 
21 namespace mfem
22 {
23 
24 class ParFiniteElementSpace;
25 
26 /// Operator that extracts Face degrees of freedom for NCMesh in parallel.
27 /** Objects of this type are typically created and owned by
28  ParFiniteElementSpace objects, see
29  ParFiniteElementSpace::GetFaceRestriction(). */
31 {
32 protected:
33  const FaceType type;
35  mutable Vector x_interp;
36 
37 public:
38  /** @brief Constructs an ParNCH1FaceRestriction.
39 
40  @param[in] fes The ParFiniteElementSpace on which this operates
41  @param[in] f_ordering Request a specific face dof ordering
42  @param[in] type Request internal or boundary faces dofs */
44  ElementDofOrdering f_ordering,
45  FaceType type);
46 
47  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
48  face E-Vector.
49 
50  @param[in] x The L-vector degrees of freedom.
51  @param[out] y The face E-Vector degrees of freedom with the given format:
52  face_dofs x vdim x nf
53  where nf is the number of interior or boundary faces
54  requested by @a type in the constructor.
55  The face_dofs are ordered according to the given
56  ElementDofOrdering. */
57  void Mult(const Vector &x, Vector &y) const override;
58 
59  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
60  L-Vector.
61 
62  @param[in] x The face E-Vector degrees of freedom with the given format:
63  face_dofs x vdim x nf
64  where nf is the number of interior or boundary faces
65  requested by @a type in the constructor.
66  The face_dofs should be ordered according to the given
67  ElementDofOrdering.
68  @param[in,out] y The L-vector degrees of freedom.
69  @param[in] a Scalar coefficient for addition. */
70  void AddMultTranspose(const Vector &x, Vector &y,
71  const double a = 1.0) const override;
72 
73  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
74  L-Vector.
75 
76  @param[in,out] x The face E-Vector degrees of freedom with the given format:
77  face_dofs x vdim x nf
78  where nf is the number of interior or boundary faces
79  requested by @a type in the constructor.
80  The face_dofs should be ordered according to the given
81  ElementDofOrdering.
82  @param[in,out] y The L-vector degrees of freedom.
83 
84  @note This method is an optimization of AddMultTranspose where the @a x
85  Vector is used and modified to avoid memory allocation and memcpy. */
86  void AddMultTransposeInPlace(Vector &x, Vector &y) const override;
87 
88 private:
89  /** @brief Compute the scatter indices: L-vector to E-vector, the offsets
90  for the gathering: E-vector to L-vector, and the interpolators from
91  coarse to fine face for master non-comforming faces.
92 
93  @param[in] f_ordering Request a specific face dof ordering.
94  @param[in] type Request internal or boundary faces dofs.
95  */
96  void ComputeScatterIndicesAndOffsets(const ElementDofOrdering f_ordering,
97  const FaceType type);
98 
99  /** @brief Compute the gather indices: E-vector to L-vector.
100 
101  Note: Requires the gather offsets to be computed.
102 
103  @param[in] f_ordering Request a specific face dof ordering.
104  @param[in] type Request internal or boundary faces dofs.
105  */
106  void ComputeGatherIndices(const ElementDofOrdering f_ordering,
107  const FaceType type);
108 
109 public: // For nvcc
110  /** @brief Apply a change of basis from coarse element basis to fine element
111  basis for the coarse face dofs.
112 
113  @param[in,out] x The dofs vector that needs coarse dofs to be express in
114  term of the fine basis.
115  */
116  void NonconformingInterpolation(Vector& x) const;
117 
118  /** @brief Apply a change of basis from fine element basis to coarse element
119  basis for the coarse face dofs.
120 
121  @param[in] x The dofs vector that needs coarse dofs to be express in term
122  of the coarse basis, the result is stored in x_interp.
123  */
124  void NonconformingTransposeInterpolation(const Vector& x) const;
125 
126  /** @brief Apply a change of basis from fine element basis to coarse element
127  basis for the coarse face dofs.
128 
129  @param[in] x The dofs vector that needs coarse dofs to be express in term
130  of the coarse basis, the result is stored in x_interp.
131  */
133 };
134 
135 /// Operator that extracts Face degrees of freedom in parallel.
136 /** Objects of this type are typically created and owned by
137  ParFiniteElementSpace objects, see
138  ParFiniteElementSpace::GetFaceRestriction(). */
140 {
141 protected:
142  /** @brief Constructs an ParL2FaceRestriction.
143 
144  @param[in] fes The ParFiniteElementSpace on which this operates
145  @param[in] f_ordering Request a specific face dof ordering
146  @param[in] type Request internal or boundary faces dofs
147  @param[in] m Request the face dofs for elem1, or both elem1 and
148  elem2
149  @param[in] build Request the ParL2FaceRestriction to compute the
150  scatter/gather indices. False should only be used
151  when inheriting from ParL2FaceRestriction. */
153  ElementDofOrdering f_ordering,
154  FaceType type,
155  L2FaceValues m,
156  bool build);
157 
158 public:
159  /** @brief Constructs an ParL2FaceRestriction.
160 
161  @param[in] fes The ParFiniteElementSpace on which this operates
162  @param[in] f_ordering Request a specific face dof ordering
163  @param[in] type Request internal or boundary faces dofs
164  @param[in] m Request the face dofs for elem1, or both elem1 and
165  elem2 */
167  ElementDofOrdering f_ordering,
168  FaceType type,
170 
171  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
172  face E-Vector.
173 
174  @param[in] x The L-vector degrees of freedom.
175  @param[out] y The face E-Vector degrees of freedom with the given format:
176  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
177  if L2FacesValues::SingleValued (face_dofs x vdim x nf),
178  where nf is the number of interior or boundary faces
179  requested by @a type in the constructor.
180  The face_dofs are ordered according to the given
181  ElementDofOrdering. */
182  void Mult(const Vector &x, Vector &y) const override;
183 
184  /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
185  given by this ParL2FaceRestriction.
186 
187  @param[in,out] mat The sparse matrix for which we want to initialize the
188  row offsets.
189  @param[in] keep_nbr_block When set to true the SparseMatrix will
190  include the rows (in addition to the columns)
191  corresponding to face-neighbor dofs. The
192  default behavior is to disregard those rows. */
193  void FillI(SparseMatrix &mat,
194  const bool keep_nbr_block = false) const override;
195 
196  /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
197  given by this ParL2FaceRestriction. @a mat contains the interior dofs
198  contribution, the @a face_mat contains the shared dofs contribution.*/
199  void FillI(SparseMatrix &mat,
200  SparseMatrix &face_mat) const;
201 
202  /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
203  pattern given by this ParL2FaceRestriction, and the values of ea_data.
204  @a mat contains the interior dofs contribution, the @a face_mat contains
205  the shared dofs contribution.*/
206  void FillJAndData(const Vector &ea_data,
207  SparseMatrix &mat,
208  SparseMatrix &face_mat) const;
209 
210  /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
211  the sparsity pattern given by this ParL2FaceRestriction, and the values of
212  fea_data.
213 
214  @param[in] fea_data The dense matrices representing the local operators
215  on each face. The format is:
216  face_dofs x face_dofs x 2 x nf.
217  On each face the first local matrix corresponds to
218  the contribution of elem1 on elem2, and the second to
219  the contribution of elem2 on elem1.
220  @param[in,out] mat The sparse matrix that is getting filled.
221  @param[in] keep_nbr_block When set to true the SparseMatrix will
222  include the rows (in addition to the columns)
223  corresponding to face-neighbor dofs. The
224  default behavior is to disregard those rows. */
225  void FillJAndData(const Vector &fea_data,
226  SparseMatrix &mat,
227  const bool keep_nbr_block = false) const override;
228 
229 private:
230  /** @brief Compute the scatter indices: L-vector to E-vector, and the offsets
231  for the gathering: E-vector to L-vector.
232 
233  @param[in] f_ordering Request a specific face dof ordering.
234  @param[in] type Request internal or boundary faces dofs.
235  */
236  void ComputeScatterIndicesAndOffsets(const ElementDofOrdering f_ordering,
237  const FaceType type);
238 
239  /** @brief Compute the gather indices: E-vector to L-vector.
240 
241  Note: Requires the gather offsets to be computed.
242 
243  @param[in] f_ordering Request a specific face dof ordering.
244  @param[in] type Request internal or boundary faces dofs.
245  */
246  void ComputeGatherIndices(const ElementDofOrdering f_ordering,
247  const FaceType type);
248 
249 public:
250  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
251  face E-Vector. Should only be used with conforming faces and when:
252  m == L2FacesValues::DoubleValued
253 
254  @param[in] x The L-vector degrees of freedom.
255  @param[out] y The face E-Vector degrees of freedom with the given format:
256  face_dofs x vdim x 2 x nf
257  where nf is the number of interior or boundary faces
258  requested by @a type in the constructor.
259  The face_dofs are ordered according to the given
260  ElementDofOrdering. */
261  void DoubleValuedConformingMult(const Vector& x, Vector& y) const override;
262 };
263 
264 /// Operator that extracts Face degrees of freedom for NCMesh in parallel.
265 /** Objects of this type are typically created and owned by
266  ParFiniteElementSpace objects, see
267  ParFiniteElementSpace::GetFaceRestriction(). */
270 {
271 public:
272  /** @brief Constructs an ParNCL2FaceRestriction.
273 
274  @param[in] fes The ParFiniteElementSpace on which this operates
275  @param[in] f_ordering Request a specific ordering
276  @param[in] type Request internal or boundary faces dofs
277  @param[in] m Request the face dofs for elem1, or both elem1 and
278  elem2 */
280  ElementDofOrdering f_ordering,
281  FaceType type,
283 
284  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
285  face E-Vector.
286 
287  @param[in] x The L-vector degrees of freedom.
288  @param[out] y The face E-Vector degrees of freedom with the given format:
289  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
290  if L2FacesValues::SingleValued (face_dofs x vdim x nf),
291  where nf is the number of interior or boundary faces
292  requested by @a type in the constructor.
293  The face_dofs are ordered according to the given
294  ElementDofOrdering. */
295  void Mult(const Vector &x, Vector &y) const override;
296 
297  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
298  L-Vector.
299 
300  @param[in] x The face E-Vector degrees of freedom with the given format:
301  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
302  if L2FacesValues::SingleValued (face_dofs x vdim x nf),
303  where nf is the number of interior or boundary faces
304  requested by @a type in the constructor.
305  The face_dofs should be ordered according to the given
306  ElementDofOrdering
307  @param[in,out] y The L-vector degrees of freedom.
308  @param[in] a Scalar coefficient for addition. */
309  void AddMultTranspose(const Vector &x, Vector &y,
310  const double a = 1.0) const override;
311 
312  /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
313  L-Vector.
314 
315  @param[in,out] x The face E-Vector degrees of freedom with the given format:
316  if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
317  if L2FacesValues::SingleValued (face_dofs x vdim x nf),
318  where nf is the number of interior or boundary faces
319  requested by @a type in the constructor.
320  The face_dofs should be ordered according to the given
321  ElementDofOrdering
322  @param[in,out] y The L-vector degrees of freedom.
323 
324  @note @a x is used for computation. */
325  void AddMultTransposeInPlace(Vector &x, Vector &y) const override;
326 
327  /** @brief Fill the I array of SparseMatrix corresponding to the sparsity
328  pattern given by this ParNCL2FaceRestriction.
329 
330  @param[in,out] mat The sparse matrix for which we want to initialize the
331  row offsets.
332  @param[in] keep_nbr_block When set to true the SparseMatrix will
333  include the rows (in addition to the columns)
334  corresponding to face-neighbor dofs. The
335  default behavior is to disregard those rows.
336 
337  @warning This method is not implemented yet. */
338  void FillI(SparseMatrix &mat,
339  const bool keep_nbr_block = false) const override;
340 
341  /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
342  given by this ParNCL2FaceRestriction. @a mat contains the interior dofs
343  contribution, the @a face_mat contains the shared dofs contribution.
344 
345  @warning This method is not implemented yet. */
346  void FillI(SparseMatrix &mat,
347  SparseMatrix &face_mat) const;
348 
349  /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
350  pattern given by this ParNCL2FaceRestriction, and the values of ea_data.
351  @a mat contains the interior dofs contribution, the @a face_mat contains
352  the shared dofs contribution.
353 
354  @warning This method is not implemented yet. */
355  void FillJAndData(const Vector &fea_data,
356  SparseMatrix &mat,
357  SparseMatrix &face_mat) const;
358 
359  /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
360  the sparsity pattern given by this ParNCL2FaceRestriction, and the values
361  of ea_data.
362 
363  @param[in] fea_data The dense matrices representing the local operators
364  on each face. The format is:
365  face_dofs x face_dofs x 2 x nf.
366  On each face the first local matrix corresponds to
367  the contribution of elem1 on elem2, and the second to
368  the contribution of elem2 on elem1.
369  @param[in,out] mat The sparse matrix that is getting filled.
370  @param[in] keep_nbr_block When set to true the SparseMatrix will
371  include the rows (in addition to the columns)
372  corresponding to face-neighbor dofs. The
373  default behavior is to disregard those rows.
374 
375  @warning This method is not implemented yet. */
376  void FillJAndData(const Vector &fea_data,
377  SparseMatrix &mat,
378  const bool keep_nbr_block = false) const override;
379 
380 private:
381  /** @brief Compute the scatter indices: L-vector to E-vector, the offsets
382  for the gathering: E-vector to L-vector, and the interpolators from
383  coarse to fine face for master non-comforming faces.
384 
385  @param[in] f_ordering Request a specific face dof ordering.
386  @param[in] type Request internal or boundary faces dofs.
387  */
388  void ComputeScatterIndicesAndOffsets(const ElementDofOrdering f_ordering,
389  const FaceType type);
390 
391  /** @brief Compute the gather indices: E-vector to L-vector.
392 
393  Note: Requires the gather offsets to be computed.
394 
395  @param[in] f_ordering Request a specific face dof ordering.
396  @param[in] type Request internal or boundary faces dofs.
397  */
398  void ComputeGatherIndices(const ElementDofOrdering f_ordering,
399  const FaceType type);
400 
401 public:
402  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
403  face E-Vector. Should only be used with nonconforming faces and when:
404  L2FaceValues m == L2FaceValues::SingleValued
405 
406  @param[in] x The L-vector degrees of freedom.
407  @param[out] y The face E-Vector degrees of freedom with the given format:
408  (face_dofs x vdim x nf),
409  where nf is the number of interior or boundary faces
410  requested by @a type in the constructor.
411  The face_dofs are ordered according to the given
412  ElementDofOrdering. */
413  void SingleValuedNonconformingMult(const Vector& x, Vector& y) const;
414 
415  /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
416  face E-Vector. Should only be used with nonconforming faces and when:
417  L2FaceValues m == L2FaceValues::DoubleValued
418 
419  @param[in] x The L-vector degrees of freedom.
420  @param[out] y The face E-Vector degrees of freedom with the given format:
421  (face_dofs x vdim x 2 x nf),
422  where nf is the number of interior or boundary faces
423  requested by @a type in the constructor.
424  The face_dofs are ordered according to the given
425  ElementDofOrdering. */
426  void DoubleValuedNonconformingMult(const Vector& x, Vector& y) const override;
427 };
428 
429 }
430 
431 #endif // MFEM_USE_MPI
432 
433 #endif // MFEM_PRESTRICTION
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
Operator that extracts Face degrees of freedom in parallel.
void DoubleValuedNonconformingMult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with no...
void FillJAndData(const Vector &fea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
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 NonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
ParNCH1FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type)
Constructs an ParNCH1FaceRestriction.
InterpolationManager interpolations
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
const L2FaceValues m
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
ParL2FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m, bool build)
Constructs an ParL2FaceRestriction.
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void SingleValuedNonconformingMult(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 AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
FaceType
Definition: mesh.hpp:45
const FiniteElementSpace & fes
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
Data type sparse matrix.
Definition: sparsemat.hpp:50
Operator that extracts Face degrees of freedom for L2 spaces.
Operator that extracts Face degrees of freedom for NCMesh in parallel.
const FiniteElementSpace & fes
This class manages the storage and computation of the interpolations from master (coarse) face to sla...
void DoubleValuedConformingMult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
void FillJAndData(const Vector &ea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
double a
Definition: lissajous.cpp:41
Operator that extracts face degrees of freedom for H1, ND, or RT FiniteElementSpaces.
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
void NonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
void NonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs...
Operator that extracts Face degrees of freedom for NCMesh in parallel.
Vector data type.
Definition: vector.hpp:58
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
ParNCL2FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m=L2FaceValues::DoubleValued)
Constructs an ParNCL2FaceRestriction.
void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const override
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 ParNCL2FaceRestr...