MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
prestriction.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_PRESTRICTION
13#define MFEM_PRESTRICTION
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_MPI
18
19#include "restriction.hpp"
20
21namespace mfem
22{
23
24class 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{
32protected:
36
37public:
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,
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 real_t 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
88private:
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
109public: // 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{
141protected:
143
144 /** @brief Constructs an ParL2FaceRestriction.
145
146 @param[in] pfes_ The ParFiniteElementSpace on which this operates
147 @param[in] f_ordering Request a specific face dof ordering
148 @param[in] type Request internal or boundary faces dofs
149 @param[in] m Request the face dofs for elem1, or both elem1 and
150 elem2
151 @param[in] build Request the ParL2FaceRestriction to compute the
152 scatter/gather indices. False should only be used
153 when inheriting from ParL2FaceRestriction. */
155 ElementDofOrdering f_ordering,
158 bool build);
159
160public:
161 /** @brief Constructs an ParL2FaceRestriction.
162
163 @param[in] fes The ParFiniteElementSpace on which this operates
164 @param[in] f_ordering Request a specific face dof ordering
165 @param[in] type Request internal or boundary faces dofs
166 @param[in] m Request the face dofs for elem1, or both elem1 and
167 elem2 */
169 ElementDofOrdering f_ordering,
172
173 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
174 face E-Vector.
175
176 @param[in] x The L-vector degrees of freedom.
177 @param[out] y The face E-Vector degrees of freedom with the given format:
178 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
179 if L2FacesValues::SingleValued (face_dofs x vdim x nf),
180 where nf is the number of interior or boundary faces
181 requested by @a type in the constructor.
182 The face_dofs are ordered according to the given
183 ElementDofOrdering. */
184 void Mult(const Vector &x, Vector &y) const override;
185
186 /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
187 given by this ParL2FaceRestriction.
188
189 @param[in,out] mat The sparse matrix for which we want to initialize the
190 row offsets.
191 @param[in] keep_nbr_block When set to true the SparseMatrix will
192 include the rows (in addition to the columns)
193 corresponding to face-neighbor dofs. The
194 default behavior is to disregard those rows. */
195 void FillI(SparseMatrix &mat,
196 const bool keep_nbr_block = false) const override;
197
198 /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
199 given by this ParL2FaceRestriction. @a mat contains the interior dofs
200 contribution, the @a face_mat contains the shared dofs contribution.*/
201 void FillI(SparseMatrix &mat,
202 SparseMatrix &face_mat) const;
203
204 /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
205 pattern given by this ParL2FaceRestriction, and the values of ea_data.
206 @a mat contains the interior dofs contribution, the @a face_mat contains
207 the shared dofs contribution.*/
208 void FillJAndData(const Vector &ea_data,
209 SparseMatrix &mat,
210 SparseMatrix &face_mat) const;
211
212 /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
213 the sparsity pattern given by this ParL2FaceRestriction, and the values of
214 fea_data.
215
216 @param[in] fea_data The dense matrices representing the local operators
217 on each face. The format is:
218 face_dofs x face_dofs x 2 x nf.
219 On each face the first local matrix corresponds to
220 the contribution of elem1 on elem2, and the second to
221 the contribution of elem2 on elem1.
222 @param[in,out] mat The sparse matrix that is getting filled.
223 @param[in] keep_nbr_block When set to true the SparseMatrix will
224 include the rows (in addition to the columns)
225 corresponding to face-neighbor dofs. The
226 default behavior is to disregard those rows. */
227 void FillJAndData(const Vector &fea_data,
228 SparseMatrix &mat,
229 const bool keep_nbr_block = false) const override;
230
231private:
232 /** @brief Compute the scatter indices: L-vector to E-vector, and the offsets
233 for the gathering: E-vector to L-vector.
234 */
235 void ComputeScatterIndicesAndOffsets();
236
237 /** @brief Compute the gather indices: E-vector to L-vector.
238
239 Note: Requires the gather offsets to be computed.
240 */
241 void ComputeGatherIndices();
242
243public:
244 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
245 face E-Vector. Should only be used with conforming faces and when:
246 m == L2FacesValues::DoubleValued
247
248 @param[in] x The L-vector degrees of freedom.
249 @param[out] y The face E-Vector degrees of freedom with the given format:
250 face_dofs x vdim x 2 x nf
251 where nf is the number of interior or boundary faces
252 requested by @a type in the constructor.
253 The face_dofs are ordered according to the given
254 ElementDofOrdering. */
255 void DoubleValuedConformingMult(const Vector& x, Vector& y) const override;
256};
257
258/// Operator that extracts Face degrees of freedom for NCMesh in parallel.
259/** Objects of this type are typically created and owned by
260 ParFiniteElementSpace objects, see
261 ParFiniteElementSpace::GetFaceRestriction(). */
264{
265public:
266 /** @brief Constructs an ParNCL2FaceRestriction.
267
268 @param[in] fes The ParFiniteElementSpace on which this operates
269 @param[in] f_ordering Request a specific ordering
270 @param[in] type Request internal or boundary faces dofs
271 @param[in] m Request the face dofs for elem1, or both elem1 and
272 elem2 */
274 ElementDofOrdering f_ordering,
277
278 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
279 face E-Vector.
280
281 @param[in] x The L-vector degrees of freedom.
282 @param[out] y The face E-Vector degrees of freedom with the given format:
283 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
284 if L2FacesValues::SingleValued (face_dofs x vdim x nf),
285 where nf is the number of interior or boundary faces
286 requested by @a type in the constructor.
287 The face_dofs are ordered according to the given
288 ElementDofOrdering. */
289 void Mult(const Vector &x, Vector &y) const override;
290
291 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
292 L-Vector.
293
294 @param[in] x The face E-Vector degrees of freedom with the given format:
295 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
296 if L2FacesValues::SingleValued (face_dofs x vdim x nf),
297 where nf is the number of interior or boundary faces
298 requested by @a type in the constructor.
299 The face_dofs should be ordered according to the given
300 ElementDofOrdering
301 @param[in,out] y The L-vector degrees of freedom.
302 @param[in] a Scalar coefficient for addition. */
303 void AddMultTranspose(const Vector &x, Vector &y,
304 const real_t a = 1.0) const override;
305
306 /** @brief Gather the degrees of freedom, i.e. goes from face E-Vector to
307 L-Vector.
308
309 @param[in,out] x The face E-Vector degrees of freedom with the given format:
310 if L2FacesValues::DoubleValued (face_dofs x vdim x 2 x nf),
311 if L2FacesValues::SingleValued (face_dofs x vdim x nf),
312 where nf is the number of interior or boundary faces
313 requested by @a type in the constructor.
314 The face_dofs should be ordered according to the given
315 ElementDofOrdering
316 @param[in,out] y The L-vector degrees of freedom.
317
318 @note @a x is used for computation. */
319 void AddMultTransposeInPlace(Vector &x, Vector &y) const override;
320
321 /** @brief Fill the I array of SparseMatrix corresponding to the sparsity
322 pattern given by this ParNCL2FaceRestriction.
323
324 @param[in,out] mat The sparse matrix for which we want to initialize the
325 row offsets.
326 @param[in] keep_nbr_block When set to true the SparseMatrix will
327 include the rows (in addition to the columns)
328 corresponding to face-neighbor dofs. The
329 default behavior is to disregard those rows.
330
331 @warning This method is not implemented yet. */
332 void FillI(SparseMatrix &mat,
333 const bool keep_nbr_block = false) const override;
334
335 /** Fill the I array of SparseMatrix corresponding to the sparsity pattern
336 given by this ParNCL2FaceRestriction. @a mat contains the interior dofs
337 contribution, the @a face_mat contains the shared dofs contribution.
338
339 @warning This method is not implemented yet. */
340 void FillI(SparseMatrix &mat,
341 SparseMatrix &face_mat) const;
342
343 /** Fill the J and Data arrays of SparseMatrix corresponding to the sparsity
344 pattern given by this ParNCL2FaceRestriction, and the values of ea_data.
345 @a mat contains the interior dofs contribution, the @a face_mat contains
346 the shared dofs contribution.
347
348 @warning This method is not implemented yet. */
349 void FillJAndData(const Vector &fea_data,
350 SparseMatrix &mat,
351 SparseMatrix &face_mat) const;
352
353 /** @brief Fill the J and Data arrays of the SparseMatrix corresponding to
354 the sparsity pattern given by this ParNCL2FaceRestriction, and the values
355 of ea_data.
356
357 @param[in] fea_data The dense matrices representing the local operators
358 on each face. The format is:
359 face_dofs x face_dofs x 2 x nf.
360 On each face the first local matrix corresponds to
361 the contribution of elem1 on elem2, and the second to
362 the contribution of elem2 on elem1.
363 @param[in,out] mat The sparse matrix that is getting filled.
364 @param[in] keep_nbr_block When set to true the SparseMatrix will
365 include the rows (in addition to the columns)
366 corresponding to face-neighbor dofs. The
367 default behavior is to disregard those rows.
368
369 @warning This method is not implemented yet. */
370 void FillJAndData(const Vector &fea_data,
371 SparseMatrix &mat,
372 const bool keep_nbr_block = false) const override;
373
374private:
375 /** @brief Compute the scatter indices: L-vector to E-vector, the offsets
376 for the gathering: E-vector to L-vector, and the interpolators from
377 coarse to fine face for master non-comforming faces.
378 */
379 void ComputeScatterIndicesAndOffsets();
380
381 /** @brief Compute the gather indices: E-vector to L-vector.
382
383 Note: Requires the gather offsets to be computed.
384 */
385 void ComputeGatherIndices();
386
387public:
388 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
389 face E-Vector. Should only be used with nonconforming faces and when:
390 L2FaceValues m == L2FaceValues::SingleValued
391
392 @param[in] x The L-vector degrees of freedom.
393 @param[out] y The face E-Vector degrees of freedom with the given format:
394 (face_dofs x vdim x nf),
395 where nf is the number of interior or boundary faces
396 requested by @a type in the constructor.
397 The face_dofs are ordered according to the given
398 ElementDofOrdering. */
399 void SingleValuedNonconformingMult(const Vector& x, Vector& y) const;
400
401 /** @brief Scatter the degrees of freedom, i.e. goes from L-Vector to
402 face E-Vector. Should only be used with nonconforming faces and when:
403 L2FaceValues m == L2FaceValues::DoubleValued
404
405 @param[in] x The L-vector degrees of freedom.
406 @param[out] y The face E-Vector degrees of freedom with the given format:
407 (face_dofs x vdim x 2 x nf),
408 where nf is the number of interior or boundary faces
409 requested by @a type in the constructor.
410 The face_dofs are ordered according to the given
411 ElementDofOrdering. */
412 void DoubleValuedNonconformingMult(const Vector& x, Vector& y) const override;
413};
414
415}
416
417#endif // MFEM_USE_MPI
418
419#endif // MFEM_PRESTRICTION
Operator that extracts face degrees of freedom for H1, ND, or RT FiniteElementSpaces.
const FiniteElementSpace & fes
This class manages the storage and computation of the interpolations from master (coarse) face to sla...
Operator that extracts Face degrees of freedom for L2 spaces.
const L2FaceValues m
const FiniteElementSpace & fes
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
Abstract parallel finite element space.
Definition pfespace.hpp:29
Operator that extracts Face degrees of freedom in parallel.
const ParFiniteElementSpace & pfes
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 Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const override
ParL2FaceRestriction(const ParFiniteElementSpace &pfes_, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m, bool build)
Constructs an ParL2FaceRestriction.
Operator that extracts Face degrees of freedom for NCMesh in parallel.
void NonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs.
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
ParNCH1FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type)
Constructs an ParNCH1FaceRestriction.
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 Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void NonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs.
InterpolationManager interpolations
void NonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs.
Operator that extracts Face degrees of freedom for NCMesh in parallel.
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 Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
ParNCL2FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m=L2FaceValues::DoubleValued)
Constructs an ParNCL2FaceRestriction.
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 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 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...
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void FillJAndData(const Vector &fea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
real_t a
Definition lissajous.cpp:41
float real_t
Definition config.hpp:43
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:75
FaceType
Definition mesh.hpp:47