MFEM  v4.5.2
Finite element discretization library
transfer.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_TRANSFER_HPP
13 #define MFEM_TRANSFER_HPP
14 
15 #include "../linalg/linalg.hpp"
16 #include "fespace.hpp"
17 
18 #ifdef MFEM_USE_MPI
19 #include "pfespace.hpp"
20 #endif
21 
22 namespace mfem
23 {
24 
25 /** @brief Base class for transfer algorithms that construct transfer Operator%s
26  between two finite element (FE) spaces. */
27 /** Generally, the two FE spaces (domain and range) can be defined on different
28  meshes. */
30 {
31 protected:
32  FiniteElementSpace &dom_fes; ///< Domain FE space
33  FiniteElementSpace &ran_fes; ///< Range FE space
34 
35  /** @brief Desired Operator::Type for the construction of all operators
36  defined by the underlying transfer algorithm. It can be ignored by
37  derived classes. */
39 
40  OperatorHandle fw_t_oper; ///< Forward true-dof operator
41  OperatorHandle bw_t_oper; ///< Backward true-dof operator
42 
43 #ifdef MFEM_USE_MPI
44  bool parallel;
45 #endif
46  bool Parallel() const
47  {
48 #ifndef MFEM_USE_MPI
49  return false;
50 #else
51  return parallel;
52 #endif
53  }
54 
56  FiniteElementSpace &fes_out,
57  const Operator &oper,
58  OperatorHandle &t_oper);
59 
60 public:
61  /** Construct a transfer algorithm between the domain, @a dom_fes_, and
62  range, @a ran_fes_, FE spaces. */
64 
65  /// Virtual destructor
66  virtual ~GridTransfer() { }
67 
68  /** @brief Set the desired Operator::Type for the construction of all
69  operators defined by the underlying transfer algorithm. */
70  /** The default value is Operator::ANY_TYPE which typically corresponds to a
71  matrix-free operator representation. Note that derived classes are not
72  required to support this setting and can ignore it. */
73  void SetOperatorType(Operator::Type type) { oper_type = type; }
74 
75  /** @brief Return an Operator that transfers GridFunction%s from the domain
76  FE space to GridFunction%s in the range FE space. */
77  virtual const Operator &ForwardOperator() = 0;
78 
79  /** @brief Return an Operator that transfers GridFunction%s from the range FE
80  space back to GridFunction%s in the domain FE space. */
81  virtual const Operator &BackwardOperator() = 0;
82 
83  /** @brief Return an Operator that transfers true-dof Vector%s from the
84  domain FE space to true-dof Vector%s in the range FE space. */
85  /** This method is implemented in the base class, based on ForwardOperator(),
86  however, derived classes can overload the construction, if necessary. */
87  virtual const Operator &TrueForwardOperator()
88  {
90  }
91 
92  /** @brief Return an Operator that transfers true-dof Vector%s from the range
93  FE space back to true-dof Vector%s in the domain FE space. */
94  /** This method is implemented in the base class, based on
95  BackwardOperator(), however, derived classes can overload the
96  construction, if necessary. */
97  virtual const Operator &TrueBackwardOperator()
98  {
100  }
101 
102  virtual bool SupportsBackwardsOperator() const { return true; }
103 };
104 
105 
106 /** @brief Transfer data between a coarse mesh and an embedded refined mesh
107  using interpolation. */
108 /** The forward, coarse-to-fine, transfer uses nodal interpolation. The
109  backward, fine-to-coarse, transfer is defined locally (on a coarse element)
110  as B = (F^t M_f F)^{-1} F^t M_f, where F is the forward transfer matrix, and
111  M_f is a mass matrix on the union of all fine elements comprising the coarse
112  element. Note that the backward transfer operator, B, is a left inverse of
113  the forward transfer operator, F, i.e. B F = I. Both F and B are defined in
114  reference space and do not depend on the actual physical shape of the mesh
115  elements.
116 
117  It is assumed that both the coarse and the fine FiniteElementSpace%s use
118  compatible types of elements, e.g. finite elements with the same map-type
119  (VALUE, INTEGRAL, H_DIV, H_CURL - see class FiniteElement). Generally, the
120  FE spaces can have different orders, however, in order for the backward
121  operator to be well-defined, the (local) number of the fine dofs should not
122  be smaller than the number of coarse dofs. */
124 {
125 protected:
126  BilinearFormIntegrator *mass_integ; ///< Ownership depends on #own_mass_integ
127  bool own_mass_integ; ///< Ownership flag for #mass_integ
128 
129  OperatorHandle F; ///< Forward, coarse-to-fine, operator
130  OperatorHandle B; ///< Backward, fine-to-coarse, operator
131 
132 public:
134  FiniteElementSpace &fine_fes)
135  : GridTransfer(coarse_fes, fine_fes),
136  mass_integ(NULL), own_mass_integ(false)
137  { }
138 
139  virtual ~InterpolationGridTransfer();
140 
141  /** @brief Assign a mass integrator to be used in the construction of the
142  backward, fine-to-coarse, transfer operator. */
143  void SetMassIntegrator(BilinearFormIntegrator *mass_integ_,
144  bool own_mass_integ_ = true);
145 
146  virtual const Operator &ForwardOperator();
147 
148  virtual const Operator &BackwardOperator();
149 };
150 
151 
152 /** @brief Transfer data in L2 and H1 finite element spaces between a coarse
153  mesh and an embedded refined mesh using L2 projection. */
154 /** The forward, coarse-to-fine, transfer uses L2 projection. The backward,
155  fine-to-coarse, transfer is defined as B = (F^t M_f F)^{-1} F^t M_f, where F
156  is the forward transfer matrix, and M_f is the mass matrix on the coarse
157  element. For L2 spaces, M_f is the mass matrix on the union of all fine
158  elements comprising the coarse element. For H1 spaces, M_f is a diagonal
159  (lumped) mass matrix computed through row-summation. Note that the backward
160  transfer operator, B, is a left inverse of the forward transfer operator, F,
161  i.e. B F = I. Both F and B are defined in physical space and, generally for
162  L2 spaces, vary between different mesh elements.
163 
164  This class supports H1 and L2 finite element spaces. Fine meshes are a
165  uniform refinement of the coarse mesh, usually created through
166  Mesh::MakeRefined. Generally, the coarse and fine FE spaces can have
167  different orders, however, in order for the backward operator to be
168  well-defined, the number of fine dofs (in a coarse element) should not be
169  smaller than the number of coarse dofs. */
171 {
172 protected:
173  /** Abstract class representing projection operator between a high-order
174  finite element space on a coarse mesh, and a low-order finite element
175  space on a refined mesh (LOR). We assume that the low-order space,
176  fes_lor, lives on a mesh obtained by refining the mesh of the high-order
177  space, fes_ho. */
178  class L2Projection : public Operator
179  {
180  public:
181  virtual void Prolongate(const Vector& x, Vector& y) const = 0;
182  virtual void ProlongateTranspose(const Vector& x, Vector& y) const = 0;
183  /// Sets relative tolerance and absolute tolerance in preconditioned
184  /// conjugate gradient solver. Only used for H1 spaces.
185  virtual void SetRelTol(double p_rtol_) = 0;
186  virtual void SetAbsTol(double p_atol_) = 0;
187  protected:
190 
192 
193  L2Projection(const FiniteElementSpace& fes_ho_,
194  const FiniteElementSpace& fes_lor_);
195 
196  void BuildHo2Lor(int nel_ho, int nel_lor,
197  const CoarseFineTransformations& cf_tr);
198 
199  void ElemMixedMass(Geometry::Type geom, const FiniteElement& fe_ho,
200  const FiniteElement& fe_lor, ElementTransformation* el_tr,
202  DenseMatrix& M_mixed_el) const;
203  };
204 
205  /** Class for projection operator between a L2 high-order finite element
206  space on a coarse mesh, and a L2 low-order finite element space on a
207  refined mesh (LOR). */
209  {
210  // The restriction and prolongation operators are represented as dense
211  // elementwise matrices (of potentially different sizes, because of mixed
212  // meshes or p-refinement). The matrix entries are stored in the R and P
213  // arrays. The entries of the i'th high-order element are stored at the
214  // index given by offsets[i].
215  mutable Array<double> R, P;
216  Array<int> offsets;
217 
218  public:
220  const FiniteElementSpace& fes_lor_);
221  /// Maps <tt>x</tt>, primal field coefficients defined on a coarse mesh
222  /// with a higher order L2 finite element space, to <tt>y</tt>, primal
223  /// field coefficients defined on a refined mesh with a low order L2
224  /// finite element space. Refined mesh should be a uniform refinement of
225  /// the coarse mesh. Coefficients are computed through minimization of L2
226  /// error between the fields.
227  virtual void Mult(const Vector& x, Vector& y) const;
228  /// Maps <tt>x</tt>, dual field coefficients defined on a refined mesh
229  /// with a low order L2 finite element space, to <tt>y</tt>, dual field
230  /// coefficients defined on a coarse mesh with a higher order L2 finite
231  /// element space. Refined mesh should be a uniform refinement of the
232  /// coarse mesh. Coefficients are computed through minimization of L2
233  /// error between the primal fields. Note, if the <tt>x</tt>-coefficients
234  /// come from ProlongateTranspose, then mass is conserved.
235  virtual void MultTranspose(const Vector& x, Vector& y) const;
236  /// Maps <tt>x</tt>, primal field coefficients defined on a refined mesh
237  /// with a low order L2 finite element space, to <tt>y</tt>, primal field
238  /// coefficients defined on a coarse mesh with a higher order L2 finite
239  /// element space. Refined mesh should be a uniform refinement of the
240  /// coarse mesh. Coefficients are computed from the mass conservative
241  /// left-inverse prolongation operation. This functionality is also
242  /// provided as an Operator by L2Prolongation.
243  virtual void Prolongate(const Vector& x, Vector& y) const;
244  /// Maps <tt>x</tt>, dual field coefficients defined on a coarse mesh with
245  /// a higher order L2 finite element space, to <tt>y</tt>, dual field
246  /// coefficients defined on a refined mesh with a low order L2 finite
247  /// element space. Refined mesh should be a uniform refinement of the
248  /// coarse mesh. Coefficients are computed from the transpose of the mass
249  /// conservative left-inverse prolongation operation. This functionality
250  /// is also provided as an Operator by L2Prolongation.
251  virtual void ProlongateTranspose(const Vector& x, Vector& y) const;
252  virtual void SetRelTol(double p_rtol_) {}
253  virtual void SetAbsTol(double p_atol_) {}
254  };
255 
256  /** Class for projection operator between a H1 high-order finite element
257  space on a coarse mesh, and a H1 low-order finite element space on a
258  refined mesh (LOR). */
260  {
261  // The restriction operator is represented as a SparseMatrix R. The
262  // prolongation operator is a dense matrix computed as the inverse of (R^T
263  // M_L R), and hence, is not stored.
264  SparseMatrix R;
265  // Used to compute P = (RTxM_LH)^(-1) M_LH^T
266  SparseMatrix M_LH;
267  SparseMatrix* RTxM_LH;
268  CGSolver pcg;
269  DSmoother Ds;
270 
271  public:
273  const FiniteElementSpace& fes_lor_);
274  virtual ~L2ProjectionH1Space();
275  /// Maps <tt>x</tt>, primal field coefficients defined on a coarse mesh
276  /// with a higher order H1 finite element space, to <tt>y</tt>, primal
277  /// field coefficients defined on a refined mesh with a low order H1
278  /// finite element space. Refined mesh should be a uniform refinement of
279  /// the coarse mesh. Coefficients are computed through minimization of L2
280  /// error between the fields.
281  virtual void Mult(const Vector& x, Vector& y) const;
282  /// Maps <tt>x</tt>, dual field coefficients defined on a refined mesh
283  /// with a low order H1 finite element space, to <tt>y</tt>, dual field
284  /// coefficients defined on a coarse mesh with a higher order H1 finite
285  /// element space. Refined mesh should be a uniform refinement of the
286  /// coarse mesh. Coefficients are computed through minimization of L2
287  /// error between the primal fields. Note, if the <tt>x</tt>-coefficients
288  /// come from ProlongateTranspose, then mass is conserved.
289  virtual void MultTranspose(const Vector& x, Vector& y) const;
290  /// Maps <tt>x</tt>, primal field coefficients defined on a refined mesh
291  /// with a low order H1 finite element space, to <tt>y</tt>, primal field
292  /// coefficients defined on a coarse mesh with a higher order H1 finite
293  /// element space. Refined mesh should be a uniform refinement of the
294  /// coarse mesh. Coefficients are computed from the mass conservative
295  /// left-inverse prolongation operation. This functionality is also
296  /// provided as an Operator by L2Prolongation.
297  virtual void Prolongate(const Vector& x, Vector& y) const;
298  /// Maps <tt>x</tt>, dual field coefficients defined on a coarse mesh with
299  /// a higher order H1 finite element space, to <tt>y</tt>, dual field
300  /// coefficients defined on a refined mesh with a low order H1 finite
301  /// element space. Refined mesh should be a uniform refinement of the
302  /// coarse mesh. Coefficients are computed from the transpose of the mass
303  /// conservative left-inverse prolongation operation. This functionality
304  /// is also provided as an Operator by L2Prolongation.
305  virtual void ProlongateTranspose(const Vector& x, Vector& y) const;
306  virtual void SetRelTol(double p_rtol_);
307  virtual void SetAbsTol(double p_atol_);
308  private:
309  /// Computes sparsity pattern and initializes R matrix. Based on
310  /// BilinearForm::AllocMat() except maps between HO elements and LOR
311  /// elements.
312  void AllocR();
313  };
314 
315  /** Mass-conservative prolongation operator going in the opposite direction
316  as L2Projection. This operator is a left inverse to the L2Projection. */
317  class L2Prolongation : public Operator
318  {
319  const L2Projection &l2proj;
320 
321  public:
322  L2Prolongation(const L2Projection &l2proj_)
323  : Operator(l2proj_.Width(), l2proj_.Height()), l2proj(l2proj_) { }
324  void Mult(const Vector &x, Vector &y) const
325  {
326  l2proj.Prolongate(x, y);
327  }
328  void MultTranspose(const Vector &x, Vector &y) const
329  {
330  l2proj.ProlongateTranspose(x, y);
331  }
332  virtual ~L2Prolongation() { }
333  };
334 
335  L2Projection *F; ///< Forward, coarse-to-fine, operator
336  L2Prolongation *B; ///< Backward, fine-to-coarse, operator
338 
339 public:
341  FiniteElementSpace &fine_fes_,
342  bool force_l2_space_ = false)
343  : GridTransfer(coarse_fes_, fine_fes_),
344  F(NULL), B(NULL), force_l2_space(force_l2_space_)
345  { }
346  virtual ~L2ProjectionGridTransfer();
347 
348  virtual const Operator &ForwardOperator();
349 
350  virtual const Operator &BackwardOperator();
351 
352  virtual bool SupportsBackwardsOperator() const;
353 private:
354  void BuildF();
355 };
356 
357 /// Matrix-free transfer operator between finite element spaces
359 {
360 private:
361  Operator* opr;
362 
363 public:
364  /// Constructs a transfer operator from \p lFESpace to \p hFESpace.
365  /** No matrices are assembled, only the action to a vector is being computed.
366  If both spaces' FE collection pointers are pointing to the same
367  collection we assume that the grid was refined while keeping the order
368  constant. If the FE collections are different, it is assumed that both
369  spaces have are using the same mesh. If the first element of the
370  high-order space is a `TensorBasisElement`, the optimized tensor-product
371  transfers are used. If not, the general transfers used. */
372  TransferOperator(const FiniteElementSpace& lFESpace,
373  const FiniteElementSpace& hFESpace);
374 
375  /// Destructor
376  virtual ~TransferOperator();
377 
378  /// @brief Interpolation or prolongation of a vector \p x corresponding to
379  /// the coarse space to the vector \p y corresponding to the fine space.
380  virtual void Mult(const Vector& x, Vector& y) const override;
381 
382  /// Restriction by applying the transpose of the Mult method.
383  /** The vector \p x corresponding to the fine space is restricted to the
384  vector \p y corresponding to the coarse space. */
385  virtual void MultTranspose(const Vector& x, Vector& y) const override;
386 };
387 
388 /// Matrix-free transfer operator between finite element spaces on the same mesh
390 {
391 private:
392  const FiniteElementSpace& lFESpace;
393  const FiniteElementSpace& hFESpace;
394  bool isvar_order;
395 
396 public:
397  /// @brief Constructs a transfer operator from \p lFESpace to \p hFESpace
398  /// which have different FE collections.
399  /** No matrices are assembled, only the action to a vector is being computed.
400  The underlying finite elements need to implement the GetTransferMatrix
401  methods. */
403  const FiniteElementSpace& hFESpace_);
404 
405  /// Destructor
407 
408  /// @brief Interpolation or prolongation of a vector \p x corresponding to
409  /// the coarse space to the vector \p y corresponding to the fine space.
410  virtual void Mult(const Vector& x, Vector& y) const override;
411 
412  /// Restriction by applying the transpose of the Mult method.
413  /** The vector \p x corresponding to the fine space is restricted to the
414  vector \p y corresponding to the coarse space. */
415  virtual void MultTranspose(const Vector& x, Vector& y) const override;
416 };
417 
418 /// @brief Matrix-free transfer operator between finite element spaces on the
419 /// same mesh exploiting the tensor product structure of the finite elements
421 {
422 private:
423  const FiniteElementSpace& lFESpace;
424  const FiniteElementSpace& hFESpace;
425  int dim;
426  int NE;
427  int D1D;
428  int Q1D;
429  Array<double> B;
430  Array<double> Bt;
431  const Operator* elem_restrict_lex_l;
432  const Operator* elem_restrict_lex_h;
433  Vector mask;
434  mutable Vector localL;
435  mutable Vector localH;
436 
437 public:
438  /// @brief Constructs a transfer operator from \p lFESpace to \p hFESpace
439  /// which have different FE collections.
440  /** No matrices are assembled, only the action to a vector is being computed.
441  The underlying finite elements need to be of type `TensorBasisElement`. It is
442  also assumed that all the elements in the spaces are of the same type. */
444  const FiniteElementSpace& lFESpace_,
445  const FiniteElementSpace& hFESpace_);
446 
447  /// Destructor
449 
450  /// @brief Interpolation or prolongation of a vector \p x corresponding to
451  /// the coarse space to the vector \p y corresponding to the fine space.
452  virtual void Mult(const Vector& x, Vector& y) const override;
453 
454  /// Restriction by applying the transpose of the Mult method.
455  /** The vector \p x corresponding to the fine space is restricted to the
456  vector \p y corresponding to the coarse space. */
457  virtual void MultTranspose(const Vector& x, Vector& y) const override;
458 };
459 
460 /// @brief Matrix-free transfer operator between finite element spaces working
461 /// on true degrees of freedom
463 {
464 private:
465  const FiniteElementSpace& lFESpace;
466  const FiniteElementSpace& hFESpace;
467  const Operator * P = nullptr;
468  const SparseMatrix * R = nullptr;
469  TransferOperator* localTransferOperator;
470  mutable Vector tmpL;
471  mutable Vector tmpH;
472 
473 public:
474  /// @brief Constructs a transfer operator working on true degrees of freedom
475  /// from \p lFESpace to \p hFESpace
476  TrueTransferOperator(const FiniteElementSpace& lFESpace_,
477  const FiniteElementSpace& hFESpace_);
478 
479  /// Destructor
481 
482  /// @brief Interpolation or prolongation of a true dof vector \p x to a true
483  /// dof vector \p y.
484  /** The true dof vector \p x corresponding to the coarse space is restricted
485  to the true dof vector \p y corresponding to the fine space. */
486  virtual void Mult(const Vector& x, Vector& y) const override;
487 
488  /// Restriction by applying the transpose of the Mult method.
489  /** The true dof vector \p x corresponding to the fine space is restricted to
490  the true dof vector \p y corresponding to the coarse space. */
491  virtual void MultTranspose(const Vector& x, Vector& y) const override;
492 };
493 
494 } // namespace mfem
495 
496 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:232
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
Definition: transfer.cpp:34
virtual void Prolongate(const Vector &x, Vector &y) const
Definition: transfer.cpp:473
Conjugate gradient method.
Definition: solvers.hpp:493
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:1426
Data type for scaled Jacobi-type smoother of sparse matrix.
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
Definition: transfer.hpp:126
L2ProjectionH1Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition: transfer.cpp:545
virtual ~TensorProductPRefinementTransferOperator()
Destructor.
Definition: transfer.cpp:1396
L2Prolongation(const L2Projection &l2proj_)
Definition: transfer.hpp:322
virtual void Mult(const Vector &x, Vector &y) const
Definition: transfer.cpp:692
FiniteElementSpace & dom_fes
Domain FE space.
Definition: transfer.hpp:32
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition: transfer.hpp:389
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:1506
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual ~TransferOperator()
Destructor.
Definition: transfer.cpp:954
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
Definition: transfer.hpp:38
virtual void ProlongateTranspose(const Vector &x, Vector &y) const =0
~TrueTransferOperator()
Destructor.
Definition: transfer.cpp:1482
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
Definition: transfer.cpp:956
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: transfer.cpp:887
virtual void Prolongate(const Vector &x, Vector &y) const
Definition: transfer.cpp:734
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
Definition: transfer.cpp:510
void BuildHo2Lor(int nel_ho, int nel_lor, const CoarseFineTransformations &cf_tr)
Definition: transfer.cpp:238
virtual ~PRefinementTransferOperator()
Destructor.
Definition: transfer.cpp:975
virtual const Operator & TrueBackwardOperator()
Return an Operator that transfers true-dof Vectors from the range FE space back to true-dof Vectors i...
Definition: transfer.hpp:97
virtual ~GridTransfer()
Virtual destructor.
Definition: transfer.hpp:66
FiniteElementSpace & ran_fes
Range FE space.
Definition: transfer.hpp:33
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition: transfer.cpp:189
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:961
Data type sparse matrix.
Definition: sparsemat.hpp:50
virtual void SetRelTol(double p_rtol_)=0
void ElemMixedMass(Geometry::Type geom, const FiniteElement &fe_ho, const FiniteElement &fe_lor, ElementTransformation *el_tr, IntegrationPointTransformation &ip_tr, DenseMatrix &M_mixed_el) const
Definition: transfer.cpp:258
OperatorHandle F
Forward, coarse-to-fine, operator.
Definition: transfer.hpp:129
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
Definition: transfer.cpp:19
virtual void Prolongate(const Vector &x, Vector &y) const =0
TensorProductPRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
Definition: transfer.cpp:1102
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition: transfer.cpp:231
InterpolationGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
Definition: transfer.hpp:133
OperatorHandle bw_t_oper
Backward true-dof operator.
Definition: transfer.hpp:41
virtual bool SupportsBackwardsOperator() const
Definition: transfer.hpp:102
virtual void Mult(const Vector &x, Vector &y) const
Definition: transfer.cpp:406
void SetOperatorType(Operator::Type type)
Set the desired Operator::Type for the construction of all operators defined by the underlying transf...
Definition: transfer.hpp:73
OperatorHandle B
Backward, fine-to-coarse, operator.
Definition: transfer.hpp:130
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition: transfer.cpp:893
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
Definition: transfer.hpp:29
OperatorHandle fw_t_oper
Forward true-dof operator.
Definition: transfer.hpp:40
void SetMassIntegrator(BilinearFormIntegrator *mass_integ_, bool own_mass_integ_=true)
Assign a mass integrator to be used in the construction of the backward, fine-to-coarse, transfer operator.
Definition: transfer.cpp:146
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:283
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:35
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: transfer.hpp:328
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Definition: transfer.cpp:1028
Matrix-free transfer operator between finite element spaces.
Definition: transfer.hpp:358
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition: transfer.hpp:123
Matrix-free transfer operator between finite element spaces on the same mesh exploiting the tensor pr...
Definition: transfer.hpp:420
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
Definition: transfer.cpp:1400
L2Prolongation * B
Backward, fine-to-coarse, operator.
Definition: transfer.hpp:336
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
bool own_mass_integ
Ownership flag for mass_integ.
Definition: transfer.hpp:127
virtual void MultTranspose(const Vector &x, Vector &y) const
Definition: transfer.cpp:438
L2Projection * F
Forward, coarse-to-fine, operator.
Definition: transfer.hpp:335
virtual void SetAbsTol(double p_atol_)=0
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: transfer.cpp:155
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: transfer.hpp:324
virtual const Operator & TrueForwardOperator()
Return an Operator that transfers true-dof Vectors from the domain FE space to true-dof Vectors in th...
Definition: transfer.hpp:87
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
Definition: transfer.cpp:977
L2ProjectionGridTransfer(FiniteElementSpace &coarse_fes_, FiniteElementSpace &fine_fes_, bool force_l2_space_=false)
Definition: transfer.hpp:340
const FiniteElementSpace & fes_lor
Definition: transfer.hpp:189
PRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
Definition: transfer.cpp:967
TrueTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator working on true degrees of freedom from lFESpace to hFESpace...
Definition: transfer.cpp:1453
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a true dof vector x to a true dof vector y.
Definition: transfer.cpp:1487
Transfer data in L2 and H1 finite element spaces between a coarse mesh and an embedded refined mesh u...
Definition: transfer.hpp:170
bool Parallel() const
Definition: transfer.hpp:46
Vector data type.
Definition: vector.hpp:60
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
Definition: transfer.cpp:759
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:70
Abstract operator.
Definition: operator.hpp:24
L2ProjectionL2Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition: transfer.cpp:285
virtual bool SupportsBackwardsOperator() const
Definition: transfer.cpp:916
TransferOperator(const FiniteElementSpace &lFESpace, const FiniteElementSpace &hFESpace)
Constructs a transfer operator from lFESpace to hFESpace.
Definition: transfer.cpp:922
Matrix-free transfer operator between finite element spaces working on true degrees of freedom...
Definition: transfer.hpp:462
virtual void MultTranspose(const Vector &x, Vector &y) const
Definition: transfer.cpp:713