MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
transfer.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_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
22namespace 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{
31protected:
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
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
60public:
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. */
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. */
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. */
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{
125protected:
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
132public:
134 FiniteElementSpace &fine_fes)
135 : GridTransfer(coarse_fes, fine_fes),
136 mass_integ(NULL), own_mass_integ(false)
137 { }
138
140
141 /** @brief Assign a mass integrator to be used in the construction of the
142 backward, fine-to-coarse, transfer operator. */
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{
172protected:
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 /// @brief Sets relative tolerance in preconditioned conjugate gradient
184 /// solver.
185 ///
186 /// Only used for H1 spaces.
187 virtual void SetRelTol(real_t p_rtol_) = 0;
188 /// @brief Sets absolute tolerance in preconditioned conjugate gradient
189 /// solver.
190 ///
191 /// Only used for H1 spaces.
192 virtual void SetAbsTol(real_t p_atol_) = 0;
193 protected:
196
198
199 L2Projection(const FiniteElementSpace& fes_ho_,
200 const FiniteElementSpace& fes_lor_);
201
202 void BuildHo2Lor(int nel_ho, int nel_lor,
203 const CoarseFineTransformations& cf_tr);
204
205 void ElemMixedMass(Geometry::Type geom, const FiniteElement& fe_ho,
206 const FiniteElement& fe_lor, ElementTransformation* el_tr,
208 DenseMatrix& M_mixed_el) const;
209 };
210
211 /** Class for projection operator between a L2 high-order finite element
212 space on a coarse mesh, and a L2 low-order finite element space on a
213 refined mesh (LOR). */
215 {
216 // The restriction and prolongation operators are represented as dense
217 // elementwise matrices (of potentially different sizes, because of mixed
218 // meshes or p-refinement). The matrix entries are stored in the R and P
219 // arrays. The entries of the i'th high-order element are stored at the
220 // index given by offsets[i].
221 mutable Array<real_t> R, P;
222 Array<int> offsets;
223
224 public:
226 const FiniteElementSpace& fes_lor_);
227 /// Maps <tt>x</tt>, primal field coefficients defined on a coarse mesh
228 /// with a higher order L2 finite element space, to <tt>y</tt>, primal
229 /// field coefficients defined on a refined mesh with a low order L2
230 /// finite element space. Refined mesh should be a uniform refinement of
231 /// the coarse mesh. Coefficients are computed through minimization of L2
232 /// error between the fields.
233 virtual void Mult(const Vector& x, Vector& y) const;
234 /// Maps <tt>x</tt>, dual field coefficients defined on a refined mesh
235 /// with a low order L2 finite element space, to <tt>y</tt>, dual field
236 /// coefficients defined on a coarse mesh with a higher order L2 finite
237 /// element space. Refined mesh should be a uniform refinement of the
238 /// coarse mesh. Coefficients are computed through minimization of L2
239 /// error between the primal fields. Note, if the <tt>x</tt>-coefficients
240 /// come from ProlongateTranspose, then mass is conserved.
241 virtual void MultTranspose(const Vector& x, Vector& y) const;
242 /// Maps <tt>x</tt>, primal field coefficients defined on a refined mesh
243 /// with a low order L2 finite element space, to <tt>y</tt>, primal field
244 /// coefficients defined on a coarse mesh with a higher order L2 finite
245 /// element space. Refined mesh should be a uniform refinement of the
246 /// coarse mesh. Coefficients are computed from the mass conservative
247 /// left-inverse prolongation operation. This functionality is also
248 /// provided as an Operator by L2Prolongation.
249 virtual void Prolongate(const Vector& x, Vector& y) const;
250 /// Maps <tt>x</tt>, dual field coefficients defined on a coarse mesh with
251 /// a higher order L2 finite element space, to <tt>y</tt>, dual field
252 /// coefficients defined on a refined mesh with a low order L2 finite
253 /// element space. Refined mesh should be a uniform refinement of the
254 /// coarse mesh. Coefficients are computed from the transpose of the mass
255 /// conservative left-inverse prolongation operation. This functionality
256 /// is also provided as an Operator by L2Prolongation.
257 virtual void ProlongateTranspose(const Vector& x, Vector& y) const;
258 virtual void SetRelTol(real_t p_rtol_) { } ///< No-op.
259 virtual void SetAbsTol(real_t p_atol_) { } ///< No-op.
260 };
261
262 /** Projection operator between a H1 high-order finite element space on a
263 coarse mesh, and a H1 low-order finite element space on a refined mesh
264 (LOR). */
266 {
267 public:
269 const FiniteElementSpace &fes_lor_);
270#ifdef MFEM_USE_MPI
272 const ParFiniteElementSpace &pfes_lor_);
273#endif
274 /// Maps <tt>x</tt>, primal field coefficients defined on a coarse mesh
275 /// with a higher order H1 finite element space, to <tt>y</tt>, primal
276 /// field coefficients defined on a refined mesh with a low order H1
277 /// finite element space. Refined mesh should be a uniform refinement of
278 /// the coarse mesh. Coefficients are computed through minimization of L2
279 /// error between the fields.
280 virtual void Mult(const Vector& x, Vector& y) const;
281 /// Maps <tt>x</tt>, dual field coefficients defined on a refined mesh
282 /// with a low order H1 finite element space, to <tt>y</tt>, dual field
283 /// coefficients defined on a coarse mesh with a higher order H1 finite
284 /// element space. Refined mesh should be a uniform refinement of the
285 /// coarse mesh. Coefficients are computed through minimization of L2
286 /// error between the primal fields. Note, if the <tt>x</tt>-coefficients
287 /// come from ProlongateTranspose, then mass is conserved.
288 virtual void MultTranspose(const Vector& x, Vector& y) const;
289 /// Maps <tt>x</tt>, primal field coefficients defined on a refined mesh
290 /// with a low order H1 finite element space, to <tt>y</tt>, primal field
291 /// coefficients defined on a coarse mesh with a higher order H1 finite
292 /// element space. Refined mesh should be a uniform refinement of the
293 /// coarse mesh. Coefficients are computed from the mass conservative
294 /// left-inverse prolongation operation. This functionality is also
295 /// provided as an Operator by L2Prolongation.
296 virtual void Prolongate(const Vector& x, Vector& y) const;
297 /// Maps <tt>x</tt>, dual field coefficients defined on a coarse mesh with
298 /// a higher order H1 finite element space, to <tt>y</tt>, dual field
299 /// coefficients defined on a refined mesh with a low order H1 finite
300 /// element space. Refined mesh should be a uniform refinement of the
301 /// coarse mesh. Coefficients are computed from the transpose of the mass
302 /// conservative left-inverse prolongation operation. This functionality
303 /// is also provided as an Operator by L2Prolongation.
304 virtual void ProlongateTranspose(const Vector& x, Vector& y) const;
305 virtual void SetRelTol(real_t p_rtol_);
306 virtual void SetAbsTol(real_t p_atol_);
307 protected:
308 /// Sets up the PCG solver (sets parameters, operator, and preconditioner)
309 void SetupPCG();
310 /// Computes on-rank R and M_LH matrices.
311 std::pair<std::unique_ptr<SparseMatrix>,
312 std::unique_ptr<SparseMatrix>> ComputeSparseRAndM_LH();
313 /// @brief Recovers vector of tdofs given a vector of dofs and a finite
314 /// element space
315 void GetTDofs(const FiniteElementSpace& fes, const Vector& x, Vector& X) const;
316 /// Sets dof values given a vector of tdofs and a finite element space
317 void SetFromTDofs(const FiniteElementSpace& fes,
318 const Vector& X,
319 Vector& x) const;
320 /// @brief Recovers a vector of dual field coefficients on the tdofs given
321 /// a vector of dual coefficients and a finite element space
323 const Vector& x,
324 Vector& X) const;
325 /// @brief Sets dual field coefficients given a vector of dual field
326 /// coefficients on the tdofs and a finite element space
328 const Vector& X,
329 Vector& x) const;
330 /// @brief Fills the vdofs_list array with a list of vdofs for a given
331 /// vdim and a given finite element space
332 void TDofsListByVDim(const FiniteElementSpace& fes,
333 int vdim,
334 Array<int>& vdofs_list) const;
335 /// Returns the inverse of an on-rank lumped mass matrix
336 void LumpedMassInverse(Vector& ML_inv) const;
337 /// @brief Computes sparsity pattern and initializes R matrix.
338 ///
339 /// Based on BilinearForm::AllocMat(), except maps between coarse HO
340 /// elements and refined LOR elements.
341 std::unique_ptr<SparseMatrix> AllocR();
342
344 std::unique_ptr<Solver> precon;
345 // The restriction operator is represented as an Operator R. The
346 // prolongation operator is a dense matrix computed as the inverse of (R^T
347 // M_L R), and hence, is not stored.
348 std::unique_ptr<Operator> R;
349 // Used to compute P = (RT*M_LH)^(-1) M_LH^T
350 std::unique_ptr<Operator> M_LH;
351 std::unique_ptr<Operator> RTxM_LH;
352 };
353
354 /** Mass-conservative prolongation operator going in the opposite direction
355 as L2Projection. This operator is a left inverse to the L2Projection. */
357 {
358 const L2Projection &l2proj;
359
360 public:
362 : Operator(l2proj_.Width(), l2proj_.Height()), l2proj(l2proj_) { }
363 void Mult(const Vector &x, Vector &y) const
364 {
365 l2proj.Prolongate(x, y);
366 }
367 void MultTranspose(const Vector &x, Vector &y) const
368 {
369 l2proj.ProlongateTranspose(x, y);
370 }
371 virtual ~L2Prolongation() { }
372 };
373
374 L2Projection *F; ///< Forward, coarse-to-fine, operator
375 L2Prolongation *B; ///< Backward, fine-to-coarse, operator
377
378public:
380 FiniteElementSpace &fine_fes_,
381 bool force_l2_space_ = false)
382 : GridTransfer(coarse_fes_, fine_fes_),
383 F(NULL), B(NULL), force_l2_space(force_l2_space_)
384 { }
386
387 virtual const Operator &ForwardOperator();
388
389 virtual const Operator &BackwardOperator();
390
391 virtual bool SupportsBackwardsOperator() const;
392private:
393 void BuildF();
394};
395
396/// Matrix-free transfer operator between finite element spaces
398{
399private:
400 Operator* opr;
401
402public:
403 /// Constructs a transfer operator from \p lFESpace to \p hFESpace.
404 /** No matrices are assembled, only the action to a vector is being computed.
405 If both spaces' FE collection pointers are pointing to the same
406 collection we assume that the grid was refined while keeping the order
407 constant. If the FE collections are different, it is assumed that both
408 spaces have are using the same mesh. If the first element of the
409 high-order space is a `TensorBasisElement`, the optimized tensor-product
410 transfers are used. If not, the general transfers used. */
411 TransferOperator(const FiniteElementSpace& lFESpace,
412 const FiniteElementSpace& hFESpace);
413
414 /// Destructor
415 virtual ~TransferOperator();
416
417 /// @brief Interpolation or prolongation of a vector \p x corresponding to
418 /// the coarse space to the vector \p y corresponding to the fine space.
419 virtual void Mult(const Vector& x, Vector& y) const override;
420
421 /// Restriction by applying the transpose of the Mult method.
422 /** The vector \p x corresponding to the fine space is restricted to the
423 vector \p y corresponding to the coarse space. */
424 virtual void MultTranspose(const Vector& x, Vector& y) const override;
425};
426
427/// Matrix-free transfer operator between finite element spaces on the same mesh
429{
430private:
431 const FiniteElementSpace& lFESpace;
432 const FiniteElementSpace& hFESpace;
433 bool isvar_order;
434
435public:
436 /// @brief Constructs a transfer operator from \p lFESpace to \p hFESpace
437 /// which have different FE collections.
438 /** No matrices are assembled, only the action to a vector is being computed.
439 The underlying finite elements need to implement the GetTransferMatrix
440 methods. */
442 const FiniteElementSpace& hFESpace_);
443
444 /// Destructor
446
447 /// @brief Interpolation or prolongation of a vector \p x corresponding to
448 /// the coarse space to the vector \p y corresponding to the fine space.
449 virtual void Mult(const Vector& x, Vector& y) const override;
450
451 /// Restriction by applying the transpose of the Mult method.
452 /** The vector \p x corresponding to the fine space is restricted to the
453 vector \p y corresponding to the coarse space. */
454 virtual void MultTranspose(const Vector& x, Vector& y) const override;
455};
456
457/// @brief Matrix-free transfer operator between finite element spaces on the
458/// same mesh exploiting the tensor product structure of the finite elements
460{
461private:
462 const FiniteElementSpace& lFESpace;
463 const FiniteElementSpace& hFESpace;
464 int dim;
465 int NE;
466 int D1D;
467 int Q1D;
469 Array<real_t> Bt;
470 const Operator* elem_restrict_lex_l;
471 const Operator* elem_restrict_lex_h;
472 Vector mask;
473 mutable Vector localL;
474 mutable Vector localH;
475
476public:
477 /// @brief Constructs a transfer operator from \p lFESpace to \p hFESpace
478 /// which have different FE collections.
479 /** No matrices are assembled, only the action to a vector is being computed.
480 The underlying finite elements need to be of type `TensorBasisElement`. It is
481 also assumed that all the elements in the spaces are of the same type. */
483 const FiniteElementSpace& lFESpace_,
484 const FiniteElementSpace& hFESpace_);
485
486 /// Destructor
488
489 /// @brief Interpolation or prolongation of a vector \p x corresponding to
490 /// the coarse space to the vector \p y corresponding to the fine space.
491 virtual void Mult(const Vector& x, Vector& y) const override;
492
493 /// Restriction by applying the transpose of the Mult method.
494 /** The vector \p x corresponding to the fine space is restricted to the
495 vector \p y corresponding to the coarse space. */
496 virtual void MultTranspose(const Vector& x, Vector& y) const override;
497};
498
499/// @brief Matrix-free transfer operator between finite element spaces working
500/// on true degrees of freedom
502{
503private:
504 const FiniteElementSpace& lFESpace;
505 const FiniteElementSpace& hFESpace;
506 const Operator * P = nullptr;
507 const SparseMatrix * R = nullptr;
508 TransferOperator* localTransferOperator;
509 mutable Vector tmpL;
510 mutable Vector tmpH;
511
512public:
513 /// @brief Constructs a transfer operator working on true degrees of freedom
514 /// from \p lFESpace to \p hFESpace
516 const FiniteElementSpace& hFESpace_);
517
518 /// Destructor
520
521 /// @brief Interpolation or prolongation of a true dof vector \p x to a true
522 /// dof vector \p y.
523 /** The true dof vector \p x corresponding to the coarse space is restricted
524 to the true dof vector \p y corresponding to the fine space. */
525 virtual void Mult(const Vector& x, Vector& y) const override;
526
527 /// Restriction by applying the transpose of the Mult method.
528 /** The true dof vector \p x corresponding to the fine space is restricted to
529 the true dof vector \p y corresponding to the coarse space. */
530 virtual void MultTranspose(const Vector& x, Vector& y) const override;
531};
532
533} // namespace mfem
534
535#endif
Abstract base class BilinearFormIntegrator.
Conjugate gradient method.
Definition solvers.hpp:513
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
Abstract class for all finite elements.
Definition fe_base.hpp:239
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
Definition transfer.hpp:30
bool Parallel() const
Definition transfer.hpp:46
virtual ~GridTransfer()
Virtual destructor.
Definition transfer.hpp:66
FiniteElementSpace & dom_fes
Domain FE space.
Definition transfer.hpp:32
OperatorHandle fw_t_oper
Forward true-dof operator.
Definition transfer.hpp:40
FiniteElementSpace & ran_fes
Range FE space.
Definition transfer.hpp:33
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
OperatorHandle bw_t_oper
Backward true-dof operator.
Definition transfer.hpp:41
virtual bool SupportsBackwardsOperator() const
Definition transfer.hpp:102
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
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
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
Definition transfer.hpp:38
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
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
Definition transfer.cpp:34
GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
Definition transfer.cpp:19
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition transfer.hpp:124
InterpolationGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
Definition transfer.hpp:133
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
Definition transfer.hpp:126
OperatorHandle F
Forward, coarse-to-fine, operator.
Definition transfer.hpp:129
bool own_mass_integ
Ownership flag for mass_integ.
Definition transfer.hpp:127
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition transfer.cpp:189
OperatorHandle B
Backward, fine-to-coarse, operator.
Definition transfer.hpp:130
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,...
Definition transfer.cpp:146
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 SetFromTDofsTranspose(const FiniteElementSpace &fes, const Vector &X, Vector &x) const
Sets dual field coefficients given a vector of dual field coefficients on the tdofs and a finite elem...
Definition transfer.cpp:962
virtual void Prolongate(const Vector &x, Vector &y) const
Definition transfer.cpp:701
void TDofsListByVDim(const FiniteElementSpace &fes, int vdim, Array< int > &vdofs_list) const
Fills the vdofs_list array with a list of vdofs for a given vdim and a given finite element space.
Definition transfer.cpp:976
void SetupPCG()
Sets up the PCG solver (sets parameters, operator, and preconditioner)
Definition transfer.cpp:638
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
Definition transfer.cpp:730
void GetTDofs(const FiniteElementSpace &fes, const Vector &x, Vector &X) const
Recovers vector of tdofs given a vector of dofs and a finite element space.
Definition transfer.cpp:920
void SetFromTDofs(const FiniteElementSpace &fes, const Vector &X, Vector &x) const
Sets dof values given a vector of tdofs and a finite element space.
Definition transfer.cpp:934
virtual void SetAbsTol(real_t p_atol_)
Sets absolute tolerance in preconditioned conjugate gradient solver.
Definition transfer.cpp:764
L2ProjectionH1Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition transfer.cpp:545
void GetTDofsTranspose(const FiniteElementSpace &fes, const Vector &x, Vector &X) const
Recovers a vector of dual field coefficients on the tdofs given a vector of dual coefficients and a f...
Definition transfer.cpp:948
virtual void Mult(const Vector &x, Vector &y) const
Definition transfer.cpp:651
std::unique_ptr< SparseMatrix > AllocR()
Computes sparsity pattern and initializes R matrix.
void LumpedMassInverse(Vector &ML_inv) const
Returns the inverse of an on-rank lumped mass matrix.
Definition transfer.cpp:997
std::pair< std::unique_ptr< SparseMatrix >, std::unique_ptr< SparseMatrix > > ComputeSparseRAndM_LH()
Computes on-rank R and M_LH matrices.
Definition transfer.cpp:772
virtual void MultTranspose(const Vector &x, Vector &y) const
Definition transfer.cpp:676
virtual void SetRelTol(real_t p_rtol_)
Sets relative tolerance in preconditioned conjugate gradient solver.
Definition transfer.cpp:759
virtual void SetRelTol(real_t p_rtol_)
No-op.
Definition transfer.hpp:258
L2ProjectionL2Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition transfer.cpp:285
virtual void Prolongate(const Vector &x, Vector &y) const
Definition transfer.cpp:473
virtual void MultTranspose(const Vector &x, Vector &y) const
Definition transfer.cpp:438
virtual void Mult(const Vector &x, Vector &y) const
Definition transfer.cpp:406
virtual void SetAbsTol(real_t p_atol_)
No-op.
Definition transfer.hpp:259
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 void ProlongateTranspose(const Vector &x, Vector &y) const =0
virtual void SetRelTol(real_t p_rtol_)=0
Sets relative tolerance in preconditioned conjugate gradient solver.
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
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition transfer.cpp:231
virtual void SetAbsTol(real_t p_atol_)=0
Sets absolute tolerance in preconditioned conjugate gradient solver.
virtual void Prolongate(const Vector &x, Vector &y) const =0
L2Prolongation(const L2Projection &l2proj_)
Definition transfer.hpp:361
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:367
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition transfer.hpp:363
Transfer data in L2 and H1 finite element spaces between a coarse mesh and an embedded refined mesh u...
Definition transfer.hpp:171
L2Projection * F
Forward, coarse-to-fine, operator.
Definition transfer.hpp:374
virtual bool SupportsBackwardsOperator() const
L2ProjectionGridTransfer(FiniteElementSpace &coarse_fes_, FiniteElementSpace &fine_fes_, bool force_l2_space_=false)
Definition transfer.hpp:379
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
L2Prolongation * B
Backward, fine-to-coarse, operator.
Definition transfer.hpp:375
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Abstract operator.
Definition operator.hpp:25
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
Type
Enumeration defining IDs for some classes derived from Operator.
Definition operator.hpp:284
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition transfer.hpp:429
PRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
virtual ~PRefinementTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
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...
Abstract parallel finite element space.
Definition pfespace.hpp:29
Data type sparse matrix.
Definition sparsemat.hpp:51
Matrix-free transfer operator between finite element spaces on the same mesh exploiting the tensor pr...
Definition transfer.hpp:460
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...
TensorProductPRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
virtual ~TensorProductPRefinementTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Matrix-free transfer operator between finite element spaces.
Definition transfer.hpp:398
virtual ~TransferOperator()
Destructor.
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...
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
TransferOperator(const FiniteElementSpace &lFESpace, const FiniteElementSpace &hFESpace)
Constructs a transfer operator from lFESpace to hFESpace.
Matrix-free transfer operator between finite element spaces working on true degrees of freedom.
Definition transfer.hpp:502
~TrueTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
TrueTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator working on true degrees of freedom from lFESpace to hFESpace.
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.
Vector data type.
Definition vector.hpp:80
float real_t
Definition config.hpp:43
Defines the coarse-fine transformations of all fine elements.
Definition ncmesh.hpp:72