MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
fespace.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_FESPACE
13#define MFEM_FESPACE
14
15#include "../config/config.hpp"
17#include "../mesh/mesh.hpp"
18#include "fe_coll.hpp"
19#include "doftrans.hpp"
20#include "restriction.hpp"
21#include <iostream>
22#include <unordered_map>
23
24namespace mfem
25{
26
27/** @brief The ordering method used when the number of unknowns per mesh node
28 (vector dimension) is bigger than 1. */
30{
31public:
32 /// %Ordering methods:
33 enum Type
34 {
35 byNODES, /**< loop first over the nodes (inner loop) then over the vector
36 dimension (outer loop); symbolically it can be represented
37 as: XXX...,YYY...,ZZZ... */
38 byVDIM /**< loop first over the vector dimension (inner loop) then over
39 the nodes (outer loop); symbolically it can be represented
40 as: XYZ,XYZ,XYZ,... */
41 };
42
43 template <Type Ord>
44 static inline int Map(int ndofs, int vdim, int dof, int vd);
45
46 template <Type Ord>
47 static void DofsToVDofs(int ndofs, int vdim, Array<int> &dofs);
48};
49
50/// @brief Type describing possible layouts for Q-vectors.
51/// @sa QuadratureInterpolator and FaceQuadratureInterpolator.
52enum class QVectorLayout
53{
54 byNODES, ///< NQPT x VDIM x NE (values) / NQPT x VDIM x DIM x NE (grads)
55 byVDIM ///< VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
56};
57
58template <> inline int
59Ordering::Map<Ordering::byNODES>(int ndofs, int vdim, int dof, int vd)
60{
61 MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
62 return (dof >= 0) ? dof+ndofs*vd : dof-ndofs*vd;
63}
64
65template <> inline int
66Ordering::Map<Ordering::byVDIM>(int ndofs, int vdim, int dof, int vd)
67{
68 MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
69 return (dof >= 0) ? vd+vdim*dof : -1-(vd+vdim*(-1-dof));
70}
71
72
73/// Constants describing the possible orderings of the DOFs in one element.
75{
76 /// Native ordering as defined by the FiniteElement.
77 /** This ordering can be used by tensor-product elements when the
78 interpolation from the DOFs to quadrature points does not use the
79 tensor-product structure. */
80 NATIVE,
81 /// Lexicographic ordering for tensor-product FiniteElements.
82 /** This ordering can be used only with tensor-product elements. */
84};
85
86// Forward declarations
87class NURBSExtension;
88class BilinearFormIntegrator;
89class QuadratureSpace;
90class QuadratureInterpolator;
91class FaceQuadratureInterpolator;
92
93
94/** @brief Class FiniteElementSpace - responsible for providing FEM view of the
95 mesh, mainly managing the set of degrees of freedom.
96
97 @details The term "degree of freedom", or "dof" for short, can mean
98 different things in different contexts. In MFEM we use "dof" to refer to
99 four closely related types of data; @ref edof "edofs", @ref ldof "ldofs",
100 @ref tdof "tdofs", and @ref vdof "vdofs".
101
102 @anchor edof @par Element DoF:
103 %Element dofs, sometimes referred to as @b edofs, are the expansion
104 coefficients used to build the linear combination of basis functions which
105 approximate a field within one element of the computational mesh. The
106 arrangement of the element dofs is determined by the basis function and
107 element types.
108 @par
109 %Element dofs are usually accessed one element at a time but they can be
110 concatenated together into a global vector when minimizing access time is
111 crucial. The global number of element dofs is not directly available from
112 the FiniteElementSpace. It can be determined by repeatedly calling
113 FiniteElementSpace::GetElementDofs and summing the lengths of the resulting
114 @a dofs arrays.
115
116 @anchor ldof @par Local DoF:
117 Most basis function types share many of their element dofs with neighboring
118 elements. Consequently, the global @ref edof "edof" vector suggested above
119 would contain many redundant entries. One of the primary roles of the
120 FiniteElementSpace is to collapse out these redundancies and
121 define a unique ordering of the remaining degrees of freedom. The
122 collapsed set of dofs are called @b "local dofs" or @b ldofs in
123 the MFEM parlance.
124 @par
125 The term @b local in this context refers to the local rank in a parallel
126 processing environment. MFEM can, of course, be used in sequential
127 computing environments but it is designed with parallel processing in mind
128 and this terminology reflects that design focus.
129 @par
130 When running in parallel the set of local dofs contains all of the degrees
131 of freedom associated with locally owned elements. When running in serial
132 all elements are locally owned so all element dofs are represented in the
133 set of local dofs.
134 @par
135 There are two important caveats regarding local dofs. First, some basis
136 function types, Nedelec and Raviart-Thomas are the prime examples, have an
137 orientation associated with each basis function. The relative orientations
138 of such basis functions in neighboring elements can lead to shared degrees
139 of freedom with opposite signs from the point of view of these neighboring
140 elements. MFEM typically chooses the orientation of the first such shared
141 degree of freedom that it encounters as the default orientation for the
142 corresponding local dof. When this local dof is referenced by a neighboring
143 element which happens to require the opposite orientation the local dof
144 index will be returned (by calls to functions such as
145 FiniteElementSpace::GetElementDofs) as a negative integer. In such cases
146 the actual offset into the vector of local dofs is @b -index-1 and the
147 value expected by this element should have the opposite sign to the value
148 stored in the local dof vector.
149 @par
150 The second important caveat only pertains to high order Nedelec basis
151 functions when shared triangular faces are present in the mesh. In this
152 very particular case the relative orientation of the face with respect to
153 its two neighboring elements can lead to different definitions of the
154 degrees of freedom associated with the interior of the face which cannot
155 be handled by simply flipping the signs of the corresponding values. The
156 DofTransformation class is designed to manage the necessary @b edof to
157 @b ldof transformations in this case. In the majority of cases the
158 DofTransformation is unnecessary and a NULL pointer will be returned in
159 place of a pointer to this object. See DofTransformation for more
160 information.
161
162 @anchor tdof @par True DoF:
163 As the name suggests "true dofs" or @b tdofs form the minimal set of data
164 values needed (along with mesh and basis function definitions) to uniquely
165 define a finite element discretization of a field. The number of true dofs
166 determines the size of the linear systems which typically need to be solved
167 in FEM simulations.
168 @par
169 Often the true dofs and the local dofs are identical, however, there are
170 important cases where they differ significantly. The first such case is
171 related to non-conforming meshes. On non-conforming meshes it is common
172 for degrees of freedom associated with "hanging" nodes, edges, or faces to
173 be constrained by degrees of freedom associated with another mesh entity.
174 In such cases the "hanging" degrees of freedom should not be considered
175 "true" degrees of freedom since their values cannot be independently
176 assigned. For this reason the FiniteElementSpace must process these
177 constraints and define a reduced set of "true" degrees of freedom which are
178 distinct from the local degrees of freedom.
179 @par
180 The second important distinction arises in parallel processing. When
181 distributing a linear system in parallel each degree of freedom must be
182 assigned to a particular processor, its owner. From the finite element
183 point of view it is convenient to distribute a computational mesh and
184 define an owning processor for each element. Since degrees of freedom may
185 be shared between neighboring elements they may also be shared between
186 neighboring processors. Another role of the FiniteElementSpace is to
187 identify the ownership of degrees of freedom which must be shared between
188 processors. Therefore the set of "true" degrees of freedom must also remove
189 redundant degrees of freedom which are owned by other processors.
190 @par
191 To summarize the set of true degrees of freedom are those degrees of
192 freedom needed to solve a linear system representing the partial
193 differential equation being modeled. True dofs differ from "local" dofs by
194 eliminating redundancies across processor boundaries and applying
195 the constraints needed to properly define fields on non-conforming meshes.
196
197 @anchor vdof @par Vector DoF:
198 %Vector dofs or @b vdofs are related to fields which are constructed using
199 multiple copies of the same set of basis functions. A typical example would
200 be the use of three instances of the scalar H1 basis functions to
201 approximate the x, y, and z components of a displacement vector field in
202 three dimensional space as often seen in elasticity simulations.
203 @par
204 %Vector dofs do not represent a specific index space the way the three
205 previous types of dofs do. Rather they are related to modifications of
206 these other index spaces to accomodate multiple copies of the underlying
207 function spaces.
208 @par
209 When using @b vdofs, i.e. when @b vdim != 1, the FiniteElementSpace only
210 manages a single set of degrees of freedom and then uses simple rules to
211 determine the appropriate offsets into the full index spaces. Two ordering
212 rules are supported; @b byNODES and @b byVDIM. See Ordering::Type for
213 details.
214 @par
215 Clearly the notion of a @b vdof is relevant in each of the three contexts
216 mentioned above so extra care must be taken whenever @b vdim != 1 to ensure
217 that the @b edof, @b ldof, or @b tdof is being interpreted correctly.
218 */
220{
223 friend void Mesh::Swap(Mesh &, bool);
224 friend class LORBase;
225
226protected:
227 /// The mesh that FE space lives on (not owned).
229
230 /// Associated FE collection (not owned).
232
233 /// %Vector dimension (number of unknowns per degree of freedom).
234 int vdim;
235
236 /** Type of ordering of the vector dofs when #vdim > 1.
237 - Ordering::byNODES - first nodes, then vector dimension,
238 - Ordering::byVDIM - first vector dimension, then nodes */
240
241 /// Number of degrees of freedom. Number of unknowns is #ndofs * #vdim.
242 int ndofs;
243
244 /** Polynomial order for each element. If empty, all elements are assumed
245 to be of the default order (fec->GetOrder()). */
247
249 int uni_fdof; ///< # of single face DOFs if all faces uniform; -1 otherwise
250 int *bdofs; ///< internal DOFs of elements if mixed/var-order; NULL otherwise
251
252 /** Variable order spaces only: DOF assignments for edges and faces, see
253 docs in MakeDofTable. For constant order spaces the tables are empty. */
255 Table var_face_dofs; ///< NOTE: also used for spaces with mixed faces
256
257 /** Additional data for the var_*_dofs tables: individual variant orders
258 (these are basically alternate J arrays for var_edge/face_dofs). */
260
261 // precalculated DOFs for each element, boundary element, and face
262 mutable Table *elem_dof; // owned (except in NURBS FE space)
263 mutable Table *elem_fos; // face orientations by element index
264 mutable Table *bdr_elem_dof; // owned (except in NURBS FE space)
265 mutable Table *bdr_elem_fos; // bdr face orientations by bdr element index
266 mutable Table *face_dof; // owned; in var-order space contains variant 0 DOFs
267
269
272 mutable Array<int> face_to_be; // NURBS FE space only
273
276
277 /** Matrix representing the prolongation from the global conforming dofs to
278 a set of intermediate partially conforming dofs, e.g. the dofs associated
279 with a "cut" space on a non-conforming mesh. */
280 mutable std::unique_ptr<SparseMatrix> cP;
281 /// Conforming restriction matrix such that cR.cP=I.
282 mutable std::unique_ptr<SparseMatrix> cR;
283 /// A version of the conforming restriction matrix for variable-order spaces.
284 mutable std::unique_ptr<SparseMatrix> cR_hp;
285 mutable bool cP_is_set;
286 /// Operator computing the action of the transpose of the restriction.
287 mutable std::unique_ptr<Operator> R_transpose;
288
289 /// Transformation to apply to GridFunctions after space Update().
291
292 /// The element restriction operators, see GetElementRestriction().
294 /// The face restriction operators, see GetFaceRestriction().
295 using key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues>;
296 struct key_hash
297 {
298 std::size_t operator()(const key_face& k) const
299 {
300 return std::get<0>(k)
301 + 2 * (int)std::get<1>(k)
302 + 4 * (int)std::get<2>(k)
303 + 8 * (int)std::get<3>(k);
304 }
305 };
306 using map_L2F = std::unordered_map<const key_face,FaceRestriction*,key_hash>;
307 mutable map_L2F L2F;
308
312
313 /** Update counter, incremented every time the space is constructed/updated.
314 Used by GridFunctions to check if they are up to date with the space. */
316
317 /** Mesh sequence number last seen when constructing the space. The space
318 needs updating if Mesh::GetSequence() is larger than this. */
320
321 /// True if at least one element order changed (variable-order space only).
323
324 bool relaxed_hp; // see SetRelaxedHpConformity()
325
326 void UpdateNURBS();
327
328 void Construct();
329 void Destroy();
330
333
334 void BuildElementToDofTable() const;
335 void BuildBdrElementToDofTable() const;
336 void BuildFaceToDofTable() const;
337
338 /** @brief Generates partial face_dof table for a NURBS space.
339
340 The table is only defined for exterior faces that coincide with a
341 boundary. */
342 void BuildNURBSFaceToDofTable() const;
343
344 /// Bit-mask representing a set of orders needed by an edge/face.
345 typedef std::uint64_t VarOrderBits;
346 static constexpr int MaxVarOrder = 8*sizeof(VarOrderBits) - 1;
347
348 /// Return the minimum order (least significant bit set) in the bit mask.
349 static int MinOrder(VarOrderBits bits);
350
351 /// Return element order: internal version of GetElementOrder without checks.
352 int GetElementOrderImpl(int i) const;
353
354 /** In a variable order space, calculate a bitmask of polynomial orders that
355 need to be represented on each edge and face. */
357 Array<VarOrderBits> &face_orders) const;
358
359 /** Build the table var_edge_dofs (or var_face_dofs) in a variable order
360 space; return total edge/face DOFs. */
361 int MakeDofTable(int ent_dim, const Array<int> &entity_orders,
362 Table &entity_dofs, Array<char> *var_ent_order);
363
364 /// Search row of a DOF table for a DOF set of size 'ndof', return first DOF.
365 int FindDofs(const Table &var_dof_table, int row, int ndof) const;
366
367 /** In a variable order space, return edge DOFs associated with a polynomial
368 order that has 'ndof' degrees of freedom. */
369 int FindEdgeDof(int edge, int ndof) const
370 { return FindDofs(var_edge_dofs, edge, ndof); }
371
372 /// Similar to FindEdgeDof, but used for mixed meshes too.
373 int FindFaceDof(int face, int ndof) const
374 { return FindDofs(var_face_dofs, face, ndof); }
375
376 int FirstFaceDof(int face, int variant = 0) const
377 { return uni_fdof >= 0 ? face*uni_fdof : var_face_dofs.GetRow(face)[variant];}
378
379 /// Return number of possible DOF variants for edge/face (var. order spaces).
380 int GetNVariants(int entity, int index) const;
381
382 /// Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
383 int GetEntityDofs(int entity, int index, Array<int> &dofs,
384 Geometry::Type master_geom = Geometry::INVALID,
385 int variant = 0) const;
386 /// Helper to get vertex, edge or face VDOFs (entity=0,1,2 resp.).
387 int GetEntityVDofs(int entity, int index, Array<int> &dofs,
388 Geometry::Type master_geom = Geometry::INVALID,
389 int variant = 0) const;
390
391 // Get degenerate face DOFs: see explanation in method implementation.
393 Geometry::Type master_geom, int variant) const;
394
395 int GetNumBorderDofs(Geometry::Type geom, int order) const;
396
397 /// Calculate the cP and cR matrices for a nonconforming mesh.
398 void BuildConformingInterpolation() const;
399
400 static void AddDependencies(SparseMatrix& deps, Array<int>& master_dofs,
401 Array<int>& slave_dofs, DenseMatrix& I,
402 int skipfirst = 0);
403
404 static bool DofFinalizable(int dof, const Array<bool>& finalized,
405 const SparseMatrix& deps);
406
407 void AddEdgeFaceDependencies(SparseMatrix &deps, Array<int>& master_dofs,
408 const FiniteElement *master_fe,
409 Array<int> &slave_dofs, int slave_face,
410 const DenseMatrix *pm) const;
411
412 /// Replicate 'mat' in the vector dimension, according to vdim ordering mode.
413 void MakeVDimMatrix(SparseMatrix &mat) const;
414
415 /// GridFunction interpolation operator applicable after mesh refinement.
417 {
418 const FiniteElementSpace* fespace;
420 Table* old_elem_dof; // Owned.
421 Table* old_elem_fos; // Owned.
422
423 Array<StatelessDofTransformation*> old_DoFTransArray;
424 mutable DofTransformation old_DoFTrans;
425
426 void ConstructDoFTransArray();
427
428 public:
429 /** Construct the operator based on the elem_dof table of the original
430 (coarse) space. The class takes ownership of the table. */
432 Table *old_elem_dof/*takes ownership*/,
433 Table *old_elem_fos/*takes ownership*/, int old_ndofs);
435 const FiniteElementSpace *coarse_fes);
436 virtual void Mult(const Vector &x, Vector &y) const;
437 virtual void MultTranspose(const Vector &x, Vector &y) const;
438 virtual ~RefinementOperator();
439 };
440
441 /// Derefinement operator, used by the friend class InterpolationGridTransfer.
443 {
444 const FiniteElementSpace *fine_fes; // Not owned.
446 Table *coarse_elem_dof; // Owned.
447 // Table *coarse_elem_fos; // Owned.
448 Table coarse_to_fine;
449 Array<int> coarse_to_ref_type;
450 Array<Geometry::Type> ref_type_to_geom;
451 Array<int> ref_type_to_fine_elem_offset;
452
453 public:
455 const FiniteElementSpace *c_fes,
456 BilinearFormIntegrator *mass_integ);
457 virtual void Mult(const Vector &x, Vector &y) const;
458 virtual ~DerefinementOperator();
459 };
460
461 /** This method makes the same assumptions as the method:
462 void GetLocalRefinementMatrices(
463 const FiniteElementSpace &coarse_fes, Geometry::Type geom,
464 DenseTensor &localP) const
465 which is defined below. It also assumes that the coarse fes and this have
466 the same vector dimension, vdim. */
467 SparseMatrix *RefinementMatrix_main(const int coarse_ndofs,
468 const Table &coarse_elem_dof,
469 const Table *coarse_elem_fos,
470 const DenseTensor localP[]) const;
471
473 DenseTensor &localP) const;
475 DenseTensor &localR) const;
476
477 /** Calculate explicit GridFunction interpolation matrix (after mesh
478 refinement). NOTE: consider using the RefinementOperator class instead
479 of the fully assembled matrix, which can take a lot of memory. */
480 SparseMatrix* RefinementMatrix(int old_ndofs, const Table* old_elem_dof,
481 const Table* old_elem_fos);
482
483 /// Calculate GridFunction restriction matrix after mesh derefinement.
484 SparseMatrix* DerefinementMatrix(int old_ndofs, const Table* old_elem_dof,
485 const Table* old_elem_fos);
486
487 /** @brief Return in @a localP the local refinement matrices that map
488 between fespaces after mesh refinement. */
489 /** This method assumes that this->mesh is a refinement of coarse_fes->mesh
490 and that the CoarseFineTransformations of this->mesh are set accordingly.
491 Another assumption is that the FEs of this use the same MapType as the FEs
492 of coarse_fes. Finally, it assumes that the spaces this and coarse_fes are
493 NOT variable-order spaces. */
494 void GetLocalRefinementMatrices(const FiniteElementSpace &coarse_fes,
495 Geometry::Type geom,
496 DenseTensor &localP) const;
497
498 /// Help function for constructors + Load().
501 int vdim = 1, int ordering = Ordering::byNODES);
502
503 /// Updates the internal mesh pointer. @warning @a new_mesh must be
504 /// <b>topologically identical</b> to the existing mesh. Used if the address
505 /// of the Mesh object has changed, e.g. in @a Mesh::Swap.
506 virtual void UpdateMeshPointer(Mesh *new_mesh);
507
508 /// Resize the elem_order array on mesh change.
509 void UpdateElementOrders();
510
511 /// @brief Copies the prolongation and restriction matrices from @a fes.
512 ///
513 /// Used for low order preconditioning on non-conforming meshes. If the DOFs
514 /// require a permutation, it will be supplied by non-NULL @a perm. NULL @a
515 /// perm indicates that no permutation is required.
517 const Array<int> *perm);
518
519public:
520 /** @brief Default constructor: the object is invalid until initialized using
521 the method Load(). */
523
524 /** @brief Copy constructor: deep copy all data from @a orig except the Mesh,
525 the FiniteElementCollection, and some derived data. */
526 /** If the @a mesh or @a fec pointers are NULL (default), then the new
527 FiniteElementSpace will reuse the respective pointers from @a orig. If
528 any of these pointers is not NULL, the given pointer will be used instead
529 of the one used by @a orig.
530
531 @note The objects pointed to by the @a mesh and @a fec parameters must be
532 either the same objects as the ones used by @a orig, or copies of them.
533 Otherwise, the behavior is undefined.
534
535 @note Derived data objects, such as the conforming prolongation and
536 restriction matrices, and the update operator, will not be copied, even
537 if they are created in the @a orig object. */
538 FiniteElementSpace(const FiniteElementSpace &orig, Mesh *mesh = NULL,
539 const FiniteElementCollection *fec = NULL);
540
545
546 /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
547 /** @note If the pointer @a ext is NULL, this constructor is equivalent to
548 the standard constructor with the same arguments minus the
549 NURBSExtension, @a ext. */
554
555 /// Copy assignment not supported
557
558 /// Returns the mesh
559 inline Mesh *GetMesh() const { return mesh; }
560
561 const NURBSExtension *GetNURBSext() const { return NURBSext; }
564
565 bool Conforming() const { return mesh->Conforming() && cP == NULL; }
566 bool Nonconforming() const { return mesh->Nonconforming() || cP != NULL; }
567
568 /// Sets the order of the i'th finite element.
569 /** By default, all elements are assumed to be of fec->GetOrder(). Once
570 SetElementOrder is called, the space becomes a variable order space. */
571 void SetElementOrder(int i, int p);
572
573 /// Returns the order of the i'th finite element.
574 int GetElementOrder(int i) const;
575
576 /// Return the maximum polynomial order.
578 { return IsVariableOrder() ? elem_order.Max() : fec->GetOrder(); }
579
580 /// Returns true if the space contains elements of varying polynomial orders.
581 bool IsVariableOrder() const { return elem_order.Size(); }
582
583 /// The returned SparseMatrix is owned by the FiniteElementSpace.
585
586 /// The returned SparseMatrix is owned by the FiniteElementSpace.
588
589 /** Return a version of the conforming restriction matrix for variable-order
590 spaces with complex hp interfaces, where some true DOFs are not owned by
591 any elements and need to be interpolated from higher order edge/face
592 variants (see also @a SetRelaxedHpConformity()). */
593 /// The returned SparseMatrix is owned by the FiniteElementSpace.
595
596 /// The returned Operator is owned by the FiniteElementSpace.
597 virtual const Operator *GetProlongationMatrix() const
598 { return GetConformingProlongation(); }
599
600 /// Return an operator that performs the transpose of GetRestrictionOperator
601 /** The returned operator is owned by the FiniteElementSpace.
602
603 For a serial conforming space, this returns NULL, indicating the identity
604 operator.
605
606 For a parallel conforming space, this will return a matrix-free
607 (Device)ConformingProlongationOperator.
608
609 For a non-conforming mesh this will return a TransposeOperator wrapping
610 the restriction matrix. */
612
613 /// An abstract operator that performs the same action as GetRestrictionMatrix
614 /** In some cases this is an optimized matrix-free implementation. The
615 returned operator is owned by the FiniteElementSpace. */
616 virtual const Operator *GetRestrictionOperator() const
617 { return GetConformingRestriction(); }
618
619 /// The returned SparseMatrix is owned by the FiniteElementSpace.
620 virtual const SparseMatrix *GetRestrictionMatrix() const
621 { return GetConformingRestriction(); }
622
623 /// The returned SparseMatrix is owned by the FiniteElementSpace.
625 { return GetHpConformingRestriction(); }
626
627 /// Return an Operator that converts L-vectors to E-vectors.
628 /** An L-vector is a vector of size GetVSize() which is the same size as a
629 GridFunction. An E-vector represents the element-wise discontinuous
630 version of the FE space.
631
632 The layout of the E-vector is: ND x VDIM x NE, where ND is the number of
633 degrees of freedom, VDIM is the vector dimension of the FE space, and NE
634 is the number of the mesh elements.
635
636 The parameter @a e_ordering describes how the local DOFs in each element
637 should be ordered in the E-vector, see ElementDofOrdering.
638
639 For discontinuous spaces, the element restriction corresponds to a
640 permutation of the degrees of freedom, implemented by the
641 L2ElementRestriction class.
642
643 The returned Operator is owned by the FiniteElementSpace. */
645 ElementDofOrdering e_ordering) const;
646
647 /// Return an Operator that converts L-vectors to E-vectors on each face.
648 virtual const FaceRestriction *GetFaceRestriction(
649 ElementDofOrdering f_ordering, FaceType,
651
652 /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
653 quadrature point values and/or derivatives (Q-vectors). */
654 /** An E-vector represents the element-wise discontinuous version of the FE
655 space and can be obtained, for example, from a GridFunction using the
656 Operator returned by GetElementRestriction().
657
658 All elements will use the same IntegrationRule, @a ir as the target
659 quadrature points.
660
661 @note The returned pointer is shared. A good practice, before using it,
662 is to set all its properties to their expected values, as other parts of
663 the code may also change them. That is, it's good to call
664 SetOutputLayout() and DisableTensorProducts() before interpolating. */
666 const IntegrationRule &ir) const;
667
668 /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
669 quadrature point values and/or derivatives (Q-vectors). */
670 /** An E-vector represents the element-wise discontinuous version of the FE
671 space and can be obtained, for example, from a GridFunction using the
672 Operator returned by GetElementRestriction().
673
674 The target quadrature points in the elements are described by the given
675 QuadratureSpace, @a qs.
676
677 @note The returned pointer is shared. A good practice, before using it,
678 is to set all its properties to their expected values, as other parts of
679 the code may also change them. That is, it's good to call
680 SetOutputLayout() and DisableTensorProducts() before interpolating. */
682 const QuadratureSpace &qs) const;
683
684 /** @brief Return a FaceQuadratureInterpolator that interpolates E-vectors to
685 quadrature point values and/or derivatives (Q-vectors).
686
687 @note The returned pointer is shared. A good practice, before using it,
688 is to set all its properties to their expected values, as other parts of
689 the code may also change them. That is, it's good to call
690 SetOutputLayout() and DisableTensorProducts() before interpolating. */
692 const IntegrationRule &ir, FaceType type) const;
693
694 /// Returns the polynomial degree of the i'th finite element.
695 /** NOTE: it is recommended to use GetElementOrder in new code. */
696 int GetOrder(int i) const { return GetElementOrder(i); }
697
698 /** Return the order of an edge. In a variable order space, return the order
699 of a specific variant, or -1 if there are no more variants. */
700 int GetEdgeOrder(int edge, int variant = 0) const;
701
702 /// Returns the polynomial degree of the i'th face finite element
703 int GetFaceOrder(int face, int variant = 0) const;
704
705 /// Returns vector dimension.
706 inline int GetVDim() const { return vdim; }
707
708 /// @brief Returns number of degrees of freedom.
709 /// This is the number of @ref ldof "Local Degrees of Freedom"
710 inline int GetNDofs() const { return ndofs; }
711
712 /// @brief Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
713 inline int GetVSize() const { return vdim * ndofs; }
714
715 /// @brief Return the number of vector true (conforming) dofs.
716 virtual int GetTrueVSize() const { return GetConformingVSize(); }
717
718 /// Returns the number of conforming ("true") degrees of freedom
719 /// (if the space is on a nonconforming mesh with hanging nodes).
720 int GetNConformingDofs() const;
721
722 int GetConformingVSize() const { return vdim * GetNConformingDofs(); }
723
724 /// Return the ordering method.
725 inline Ordering::Type GetOrdering() const { return ordering; }
726
727 const FiniteElementCollection *FEColl() const { return fec; }
728
729 /// Number of all scalar vertex dofs
730 int GetNVDofs() const { return nvdofs; }
731 /// Number of all scalar edge-interior dofs
732 int GetNEDofs() const { return nedofs; }
733 /// Number of all scalar face-interior dofs
734 int GetNFDofs() const { return nfdofs; }
735
736 /// Returns number of vertices in the mesh.
737 inline int GetNV() const { return mesh->GetNV(); }
738
739 /// Returns number of elements in the mesh.
740 inline int GetNE() const { return mesh->GetNE(); }
741
742 /// Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
743 /** The co-dimension 1 entities are those that have dimension 1 less than the
744 mesh dimension, e.g. for a 2D mesh, the faces are the 1D entities, i.e.
745 the edges. */
746 inline int GetNF() const { return mesh->GetNumFaces(); }
747
748 /// Returns number of boundary elements in the mesh.
749 inline int GetNBE() const { return mesh->GetNBE(); }
750
751 /// Returns the number of faces according to the requested type.
752 /** If type==Boundary returns only the "true" number of boundary faces
753 contrary to GetNBE() that returns "fake" boundary faces associated to
754 visualization for GLVis.
755 Similarly, if type==Interior, the "fake" boundary faces associated to
756 visualization are counted as interior faces. */
757 inline int GetNFbyType(FaceType type) const
758 { return mesh->GetNFbyType(type); }
759
760 /// Returns the type of element i.
761 inline int GetElementType(int i) const
762 { return mesh->GetElementType(i); }
763
764 /// Returns the vertices of element i.
765 inline void GetElementVertices(int i, Array<int> &vertices) const
766 { mesh->GetElementVertices(i, vertices); }
767
768 /// Returns the type of boundary element i.
769 inline int GetBdrElementType(int i) const
770 { return mesh->GetBdrElementType(i); }
771
772 /// Returns ElementTransformation for the @a i-th element.
775
776 /** @brief Returns the transformation defining the @a i-th element in the
777 user-defined variable @a ElTr. */
780
781 /// Returns ElementTransformation for the @a i-th boundary element.
784
785 int GetAttribute(int i) const { return mesh->GetAttribute(i); }
786
787 int GetBdrAttribute(int i) const { return mesh->GetBdrAttribute(i); }
788
789 /// @anchor getdof @name Local DoF Access Members
790 /// These member functions produce arrays of local degree of freedom
791 /// indices, see @ref ldof. If @b vdim == 1 these indices can be used to
792 /// access entries in GridFunction, LinearForm, and BilinearForm objects.
793 /// If @b vdim != 1 the corresponding @ref getvdof "Get*VDofs" methods
794 /// should be used instead or one of the @ref dof2vdof "DofToVDof" methods
795 /// could be used to produce the appropriate offsets from these local dofs.
796 ///@{
797
798 /// @brief Returns indices of degrees of freedom of element 'elem'.
799 /// The returned indices are offsets into an @ref ldof vector. See also
800 /// GetElementVDofs().
801 ///
802 /// @note In many cases the returned DofTransformation object will be NULL.
803 /// In other cases see the documentation of the DofTransformation class for
804 /// guidance on its role in performing @ref edof to @ref ldof transformations
805 /// on local vectors and matrices. At present the DofTransformation is only
806 /// needed for Nedelec basis functions of order 2 and above on 3D elements
807 /// with triangular faces.
808 ///
809 /// @note The returned object should NOT be deleted by the caller.
810 DofTransformation *GetElementDofs(int elem, Array<int> &dofs) const;
811
812 /// @brief The same as GetElementDofs(), but with a user-allocated
813 /// DofTransformation object. @a doftrans must be allocated in advance and
814 /// will be owned by the caller. The user can use the
815 /// DofTransformation::GetDofTransformation method on the returned
816 /// @a doftrans object to detect if the DofTransformation should actually be
817 /// used.
818 virtual void GetElementDofs(int elem, Array<int> &dofs,
819 DofTransformation &doftrans) const;
820
821 /// @brief Returns indices of degrees of freedom for boundary element 'bel'.
822 /// The returned indices are offsets into an @ref ldof vector. See also
823 /// GetBdrElementVDofs().
824 ///
825 /// @note In many cases the returned DofTransformation object will be NULL.
826 /// In other cases see the documentation of the DofTransformation class for
827 /// guidance on its role in performing @ref edof to @ref ldof transformations
828 /// on local vectors and matrices. At present the DofTransformation is only
829 /// needed for Nedelec basis functions of order 2 and above on 3D elements
830 /// with triangular faces.
831 ///
832 /// @note The returned object should NOT be deleted by the caller.
833 DofTransformation *GetBdrElementDofs(int bel, Array<int> &dofs) const;
834
835 /// @brief The same as GetBdrElementDofs(), but with a user-allocated
836 /// DofTransformation object. @a doftrans must be allocated in advance and
837 /// will be owned by the caller. The user can use the
838 /// DofTransformation::GetDofTransformation method on the returned
839 /// @a doftrans object to detect if the DofTransformation should actually be
840 /// used.
841 virtual void GetBdrElementDofs(int bel, Array<int> &dofs,
842 DofTransformation &doftrans) const;
843
844 /// @brief Returns the indices of the degrees of freedom for the specified
845 /// face, including the DOFs for the edges and the vertices of the face.
846 ///
847 /// In variable order spaces, multiple variants of DOFs can be returned.
848 /// See GetEdgeDofs() for more details.
849 /// @return Order of the selected variant, or -1 if there are no more
850 /// variants.
851 ///
852 /// The returned indices are offsets into an @ref ldof vector. See also
853 /// GetFaceVDofs().
854 virtual int GetFaceDofs(int face, Array<int> &dofs, int variant = 0) const;
855
856 /// @brief Returns the indices of the degrees of freedom for the specified
857 /// edge, including the DOFs for the vertices of the edge.
858 ///
859 /// In variable order spaces, multiple sets of DOFs may exist on an edge,
860 /// corresponding to the different polynomial orders of incident elements.
861 /// The 'variant' parameter is the zero-based index of the desired DOF set.
862 /// The variants are ordered from lowest polynomial degree to the highest.
863 /// @return Order of the selected variant, or -1 if there are no more
864 /// variants.
865 ///
866 /// The returned indices are offsets into an @ref ldof vector. See also
867 /// GetEdgeVDofs().
868 int GetEdgeDofs(int edge, Array<int> &dofs, int variant = 0) const;
869
870 /// @brief Returns the indices of the degrees of freedom for the specified
871 /// vertices.
872 ///
873 /// The returned indices are offsets into an @ref ldof vector. See also
874 /// GetVertexVDofs().
875 void GetVertexDofs(int i, Array<int> &dofs) const;
876
877 /// @brief Returns the indices of the degrees of freedom for the interior
878 /// of the specified element.
879 ///
880 /// Specifically this refers to degrees of freedom which are not associated
881 /// with the vertices, edges, or faces of the mesh. This method may be
882 /// useful in conjunction with schemes which process shared and non-shared
883 /// degrees of freedom differently such as static condensation.
884 ///
885 /// The returned indices are offsets into an @ref ldof vector. See also
886 /// GetElementInteriorVDofs().
887 void GetElementInteriorDofs(int i, Array<int> &dofs) const;
888
889 /// @brief Returns the number of degrees of freedom associated with the
890 /// interior of the specified element.
891 ///
892 /// See GetElementInteriorDofs() for more information or to obtain the
893 /// relevant indices.
894 int GetNumElementInteriorDofs(int i) const;
895
896 /// @brief Returns the indices of the degrees of freedom for the interior
897 /// of the specified face.
898 ///
899 /// Specifically this refers to degrees of freedom which are not associated
900 /// with the vertices, edges, or cell interiors of the mesh. This method may
901 /// be useful in conjunction with schemes which process shared and non-shared
902 /// degrees of freedom differently such as static condensation.
903 ///
904 /// The returned indices are offsets into an @ref ldof vector. See also
905 /// GetFaceInteriorVDofs().
906 void GetFaceInteriorDofs(int i, Array<int> &dofs) const;
907
908 /// @brief Returns the indices of the degrees of freedom for the interior
909 /// of the specified edge.
910 ///
911 /// The returned indices are offsets into an @ref ldof vector. See also
912 /// GetEdgeInteriorVDofs().
913 void GetEdgeInteriorDofs(int i, Array<int> &dofs) const;
914 ///@}
915
916 /** @brief Returns indices of degrees of freedom for NURBS patch index
917 @a patch. Cartesian ordering is used, for the tensor-product degrees of
918 freedom. */
919 void GetPatchDofs(int patch, Array<int> &dofs) const;
920
921 /// @anchor dof2vdof @name DoF To VDoF Conversion methods
922 /// These methods convert between local dof and local vector dof using the
923 /// appropriate relationship based on the Ordering::Type defined in this
924 /// FiniteElementSpace object.
925 ///
926 /// These methods assume the index set has a range [0, GetNDofs()) which
927 /// will be mapped to the range [0, GetVSize()). This assumption can be
928 /// changed in the forward mappings by passing a value for @a ndofs which
929 /// differs from that returned by GetNDofs().
930 ///
931 /// @note These methods, with the exception of VDofToDof(), are designed to
932 /// produce the correctly encoded values when dof entries are negative,
933 /// see @ref ldof for more on negative dof indices.
934 ///
935 /// @warning When MFEM_DEBUG is enabled at build time the forward mappings
936 /// will verify that each @a dof lies in the proper range. If MFEM_DEBUG is
937 /// disabled no range checking is performed.
938 ///@{
939
940 /// @brief Returns the indices of all of the VDofs for the specified
941 /// dimension 'vd'.
942 ///
943 /// The @a ndofs parameter can be used to indicate the total number of Dofs
944 /// associated with each component of @b vdim. If @a ndofs is -1 (the
945 /// default value), then the number of Dofs is determined by the
946 /// FiniteElementSpace::GetNDofs().
947 ///
948 /// @note This method does not resize the @a dofs array. It takes the range
949 /// of dofs [0, dofs.Size()) and converts these to @ref vdof "vdofs" and
950 /// stores the results in the @a dofs array.
951 void GetVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
952
953 /// @brief Compute the full set of @ref vdof "vdofs" corresponding to each
954 /// entry in @a dofs.
955 ///
956 /// @details Produces a set of @ref vdof "vdofs" of
957 /// length @b vdim * dofs.Size() corresponding to the entries contained in
958 /// the @a dofs array.
959 ///
960 /// The @a ndofs parameter can be used to indicate the total number of Dofs
961 /// associated with each component of @b vdim. If @a ndofs is -1 (the
962 /// default value), then the number of Dofs is <determined by the
963 /// FiniteElementSpace::GetNDofs().
964 ///
965 /// @note The @a dofs array is overwritten and resized to accomodate the
966 /// new values.
967 void DofsToVDofs(Array<int> &dofs, int ndofs = -1) const;
968
969 /// @brief Compute the set of @ref vdof "vdofs" corresponding to each entry
970 /// in @a dofs for the given vector index @a vd.
971 ///
972 /// The @a ndofs parameter can be used to indicate the total number of Dofs
973 /// associated with each component of @b vdim. If @a ndofs is -1 (the
974 /// default value), then the number of Dofs is <determined by the
975 /// FiniteElementSpace::GetNDofs().
976 ///
977 /// @note The @a dofs array is overwritten with the new values but its size
978 /// will not be altered.
979 void DofsToVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
980
981 /// @brief Compute a single @ref vdof corresponding to the index @a dof and
982 /// the vector index @a vd.
983 ///
984 /// The @a ndofs parameter can be used to indicate the total number of Dofs
985 /// associated with each component of @b vdim. If @a ndofs is -1 (the
986 /// default value), then the number of Dofs is <determined by the
987 /// FiniteElementSpace::GetNDofs().
988 int DofToVDof(int dof, int vd, int ndofs = -1) const;
989
990 /// @brief Compute the inverse of the Dof to VDof mapping for a single
991 /// index @a vdof.
992 ///
993 /// @warning This method is only intended for use with positive indices.
994 /// Passing a negative value for @a vdof will produce an invalid result.
995 int VDofToDof(int vdof) const
996 { return (ordering == Ordering::byNODES) ? (vdof%ndofs) : (vdof/vdim); }
997
998 ///@}
999
1000 /// @brief Remove the orientation information encoded into an array of dofs
1001 /// Some basis function types have a relative orientation associated with
1002 /// degrees of freedom shared between neighboring elements, see @ref ldof
1003 /// for more information. An orientation mismatch is indicated in the dof
1004 /// indices by a negative index value. This method replaces such negative
1005 /// indices with the corresponding positive offsets.
1006 ///
1007 /// @note The name of this method reflects the fact that it is most often
1008 /// applied to sets of @ref vdof "Vector Dofs" but it would work equally
1009 /// well on sets of @ref ldof "Local Dofs".
1010 static void AdjustVDofs(Array<int> &vdofs);
1011
1012 /// Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
1013 static inline int EncodeDof(int entity_base, int idx)
1014 { return (idx >= 0) ? (entity_base + idx) : (-1-(entity_base + (-1-idx))); }
1015
1016 /// Helper to return the DOF associated with a sign encoded DOF
1017 static inline int DecodeDof(int dof)
1018 { return (dof >= 0) ? dof : (-1 - dof); }
1019
1020 /// Helper to determine the DOF and sign of a sign encoded DOF
1021 static inline int DecodeDof(int dof, real_t& sign)
1022 { return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof)); }
1023
1024 /// @anchor getvdof @name Local Vector DoF Access Members
1025 /// These member functions produce arrays of local vector degree of freedom
1026 /// indices, see @ref ldof and @ref vdof. These indices can be used to
1027 /// access entries in GridFunction, LinearForm, and BilinearForm objects
1028 /// regardless of the value of @b vdim.
1029 /// @{
1030
1031 /// @brief Returns indices of degrees of freedom for the @a i'th element.
1032 /// The returned indices are offsets into an @ref ldof vector with @b vdim
1033 /// not necessarily equal to 1. The returned indices are always ordered
1034 /// byNODES, irrespective of whether the space is byNODES or byVDIM.
1035 /// See also GetElementDofs().
1036 ///
1037 /// @note In many cases the returned DofTransformation object will be NULL.
1038 /// In other cases see the documentation of the DofTransformation class for
1039 /// guidance on its role in performing @ref edof to @ref ldof transformations
1040 /// on local vectors and matrices. At present the DofTransformation is only
1041 /// needed for Nedelec basis functions of order 2 and above on 3D elements
1042 /// with triangular faces.
1043 ///
1044 /// @note The returned object should NOT be deleted by the caller.
1045 DofTransformation *GetElementVDofs(int i, Array<int> &vdofs) const;
1046
1047 /// @brief The same as GetElementVDofs(), but with a user-allocated
1048 /// DofTransformation object. @a doftrans must be allocated in advance and
1049 /// will be owned by the caller. The user can use the
1050 /// DofTransformation::GetDofTransformation method on the returned
1051 /// @a doftrans object to detect if the DofTransformation should actually be
1052 /// used.
1053 void GetElementVDofs(int i, Array<int> &vdofs,
1054 DofTransformation &doftrans) const;
1055
1056 /// @brief Returns indices of degrees of freedom for @a i'th boundary
1057 /// element.
1058 /// The returned indices are offsets into an @ref ldof vector with @b vdim
1059 /// not necessarily equal to 1. See also GetBdrElementDofs().
1060 ///
1061 /// @note In many cases the returned DofTransformation object will be NULL.
1062 /// In other cases see the documentation of the DofTransformation class for
1063 /// guidance on its role in performing @ref edof to @ref ldof transformations
1064 /// on local vectors and matrices. At present the DofTransformation is only
1065 /// needed for Nedelec basis functions of order 2 and above on 3D elements
1066 /// with triangular faces.
1067 ///
1068 /// @note The returned object should NOT be deleted by the caller.
1069 DofTransformation *GetBdrElementVDofs(int i, Array<int> &vdofs) const;
1070
1071 /// @brief The same as GetBdrElementVDofs(), but with a user-allocated
1072 /// DofTransformation object. @a doftrans must be allocated in advance and
1073 /// will be owned by the caller. The user can use the
1074 /// DofTransformation::GetDofTransformation method on the returned
1075 /// @a doftrans object to detect if the DofTransformation should actually be
1076 /// used.
1077 void GetBdrElementVDofs(int i, Array<int> &vdofs,
1078 DofTransformation &doftrans) const;
1079
1080 /// Returns indices of degrees of freedom in @a vdofs for NURBS patch @a i.
1081 void GetPatchVDofs(int i, Array<int> &vdofs) const;
1082
1083 /// @brief Returns the indices of the degrees of freedom for the specified
1084 /// face, including the DOFs for the edges and the vertices of the face.
1085 ///
1086 /// The returned indices are offsets into an @ref ldof vector with @b vdim
1087 /// not necessarily equal to 1. See GetFaceDofs() for more information.
1088 void GetFaceVDofs(int i, Array<int> &vdofs) const;
1089
1090 /// @brief Returns the indices of the degrees of freedom for the specified
1091 /// edge, including the DOFs for the vertices of the edge.
1092 ///
1093 /// The returned indices are offsets into an @ref ldof vector with @b vdim
1094 /// not necessarily equal to 1. See GetEdgeDofs() for more information.
1095 void GetEdgeVDofs(int i, Array<int> &vdofs) const;
1096
1097 /// @brief Returns the indices of the degrees of freedom for the specified
1098 /// vertices.
1099 ///
1100 /// The returned indices are offsets into an @ref ldof vector with @b vdim
1101 /// not necessarily equal to 1. See also GetVertexDofs().
1102 void GetVertexVDofs(int i, Array<int> &vdofs) const;
1103
1104 /// @brief Returns the indices of the degrees of freedom for the interior
1105 /// of the specified element.
1106 ///
1107 /// The returned indices are offsets into an @ref ldof vector with @b vdim
1108 /// not necessarily equal to 1. See GetElementInteriorDofs() for more
1109 /// information.
1110 void GetElementInteriorVDofs(int i, Array<int> &vdofs) const;
1111
1112 /// @brief Returns the indices of the degrees of freedom for the interior
1113 /// of the specified edge.
1114 ///
1115 /// The returned indices are offsets into an @ref ldof vector with @b vdim
1116 /// not necessarily equal to 1. See also GetEdgeInteriorDofs().
1117 void GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const;
1118 /// @}
1119
1120 /// (@deprecated) Use the Update() method if the space or mesh changed.
1121 MFEM_DEPRECATED void RebuildElementToDofTable();
1122
1123 /** @brief Reorder the scalar DOFs based on the element ordering.
1124
1125 The new ordering is constructed as follows: 1) loop over all elements as
1126 ordered in the Mesh; 2) for each element, assign new indices to all of
1127 its current DOFs that are still unassigned; the new indices we assign are
1128 simply the sequence `0,1,2,...`; if there are any signed DOFs their sign
1129 is preserved. */
1131
1133
1134 /** @brief Return a reference to the internal Table that stores the lists of
1135 scalar dofs, for each mesh element, as returned by GetElementDofs(). */
1136 const Table &GetElementToDofTable() const { return *elem_dof; }
1137
1138 /** @brief Return a reference to the internal Table that stores the lists of
1139 scalar dofs, for each boundary mesh element, as returned by
1140 GetBdrElementDofs(). */
1143
1144 /** @brief Return a reference to the internal Table that stores the lists of
1145 scalar dofs, for each face in the mesh, as returned by GetFaceDofs(). In
1146 this context, "face" refers to a (dim-1)-dimensional mesh entity. */
1147 /** @note In the case of a NURBS space, the rows corresponding to interior
1148 faces will be empty. */
1150 { if (!face_dof) { BuildFaceToDofTable(); } return *face_dof; }
1151
1152 /** @brief Initialize internal data that enables the use of the methods
1153 GetElementForDof() and GetLocalDofForDof(). */
1154 void BuildDofToArrays();
1155
1156 /// Return the index of the first element that contains dof @a i.
1157 /** This method can be called only after setup is performed using the method
1158 BuildDofToArrays(). */
1159 int GetElementForDof(int i) const { return dof_elem_array[i]; }
1160 /// Return the local dof index in the first element that contains dof @a i.
1161 /** This method can be called only after setup is performed using the method
1162 BuildDofToArrays(). */
1163 int GetLocalDofForDof(int i) const { return dof_ldof_array[i]; }
1164
1165 /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1166 associated with i'th element in the mesh object.
1167 Note: The method has been updated to abort instead of returning NULL for
1168 an empty partition. */
1169 virtual const FiniteElement *GetFE(int i) const;
1170
1171 /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1172 associated with i'th boundary face in the mesh object. */
1173 const FiniteElement *GetBE(int i) const;
1174
1175 /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1176 associated with i'th face in the mesh object. Faces in this case refer
1177 to the MESHDIM-1 primitive so in 2D they are segments and in 1D they are
1178 points.*/
1179 const FiniteElement *GetFaceElement(int i) const;
1180
1181 /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1182 associated with i'th edge in the mesh object. */
1183 const FiniteElement *GetEdgeElement(int i, int variant = 0) const;
1184
1185 /// Return the trace element from element 'i' to the given 'geom_type'
1186 const FiniteElement *GetTraceElement(int i, Geometry::Type geom_type) const;
1187
1188 /** @brief Mark degrees of freedom associated with boundary elements with
1189 the specified boundary attributes (marked in 'bdr_attr_is_ess').
1190 For spaces with 'vdim' > 1, the 'component' parameter can be used
1191 to restricts the marked vDOFs to the specified component. */
1192 virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
1193 Array<int> &ess_vdofs,
1194 int component = -1) const;
1195
1196 /** @brief Get a list of essential true dofs, ess_tdof_list, corresponding to the
1197 boundary attributes marked in the array bdr_attr_is_ess.
1198 For spaces with 'vdim' > 1, the 'component' parameter can be used
1199 to restricts the marked tDOFs to the specified component. */
1200 virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
1201 Array<int> &ess_tdof_list,
1202 int component = -1) const;
1203
1204 /** @brief Get a list of all boundary true dofs, @a boundary_dofs. For spaces
1205 with 'vdim' > 1, the 'component' parameter can be used to restricts the
1206 marked tDOFs to the specified component. Equivalent to
1207 FiniteElementSpace::GetEssentialTrueDofs with all boundary attributes
1208 marked as essential. */
1209 void GetBoundaryTrueDofs(Array<int> &boundary_dofs, int component = -1);
1210
1211 /// Convert a Boolean marker array to a list containing all marked indices.
1212 static void MarkerToList(const Array<int> &marker, Array<int> &list);
1213
1214 /** @brief Convert an array of indices (list) to a Boolean marker array where all
1215 indices in the list are marked with the given value and the rest are set
1216 to zero. */
1217 static void ListToMarker(const Array<int> &list, int marker_size,
1218 Array<int> &marker, int mark_val = -1);
1219
1220 /** @brief For a partially conforming FE space, convert a marker array (nonzero
1221 entries are true) on the partially conforming dofs to a marker array on
1222 the conforming dofs. A conforming dofs is marked iff at least one of its
1223 dependent dofs is marked. */
1224 void ConvertToConformingVDofs(const Array<int> &dofs, Array<int> &cdofs);
1225
1226 /** @brief For a partially conforming FE space, convert a marker array (nonzero
1227 entries are true) on the conforming dofs to a marker array on the
1228 (partially conforming) dofs. A dof is marked iff it depends on a marked
1229 conforming dofs, where dependency is defined by the ConformingRestriction
1230 matrix; in other words, a dof is marked iff it corresponds to a marked
1231 conforming dof. */
1232 void ConvertFromConformingVDofs(const Array<int> &cdofs, Array<int> &dofs);
1233
1234 /** @brief Generate the global restriction matrix from a discontinuous
1235 FE space to the continuous FE space of the same polynomial degree. */
1237
1238 /** @brief Generate the global restriction matrix from a discontinuous
1239 FE space to the piecewise constant FE space. */
1241
1242 /** @brief Construct the restriction matrix from the FE space given by
1243 (*this) to the lower degree FE space given by (*lfes) which
1244 is defined on the same mesh. */
1246
1247 /** @brief Construct and return an Operator that can be used to transfer
1248 GridFunction data from @a coarse_fes, defined on a coarse mesh, to @a
1249 this FE space, defined on a refined mesh. */
1250 /** It is assumed that the mesh of this FE space is a refinement of the mesh
1251 of @a coarse_fes and the CoarseFineTransformations returned by the method
1252 Mesh::GetRefinementTransforms() of the refined mesh are set accordingly.
1253 The Operator::Type of @a T can be set to request an Operator of the set
1254 type. Currently, only Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE
1255 (matrix-free) are supported. When Operator::ANY_TYPE is requested, the
1256 choice of the particular Operator sub-class is left to the method. This
1257 method also works in parallel because the transfer operator is local to
1258 the MPI task when the input is a synchronized ParGridFunction. */
1259 void GetTransferOperator(const FiniteElementSpace &coarse_fes,
1260 OperatorHandle &T) const;
1261
1262 /** @brief Construct and return an Operator that can be used to transfer
1263 true-dof data from @a coarse_fes, defined on a coarse mesh, to @a this FE
1264 space, defined on a refined mesh.
1265
1266 This method calls GetTransferOperator() and multiplies the result by the
1267 prolongation operator of @a coarse_fes on the right, and by the
1268 restriction operator of this FE space on the left.
1269
1270 The Operator::Type of @a T can be set to request an Operator of the set
1271 type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and
1272 Operator::ANY_TYPE (matrix-free). In parallel, the supported types are:
1273 Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated
1274 as Operator::ANY_TYPE: the operator representation choice is made by this
1275 method. */
1276 virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
1277 OperatorHandle &T) const;
1278
1279 /** @brief Reflect changes in the mesh: update number of DOFs, etc. Also, calculate
1280 GridFunction transformation operator (unless want_transform is false).
1281 Safe to call multiple times, does nothing if space already up to date. */
1282 virtual void Update(bool want_transform = true);
1283
1284 /// Get the GridFunction update operator.
1285 const Operator* GetUpdateOperator() { Update(); return Th.Ptr(); }
1286
1287 /// Return the update operator in the given OperatorHandle, @a T.
1289
1290 /** @brief Set the ownership of the update operator: if set to false, the
1291 Operator returned by GetUpdateOperator() must be deleted outside the
1292 FiniteElementSpace. */
1293 /** The update operator ownership is automatically reset to true when a new
1294 update operator is created by the Update() method. */
1296
1297 /// Specify the Operator::Type to be used by the update operators.
1298 /** The default type is Operator::ANY_TYPE which leaves the choice to this
1299 class. The other currently supported option is Operator::MFEM_SPARSEMAT
1300 which is only guaranteed to be honored for a refinement update operator.
1301 Any other type will be treated as Operator::ANY_TYPE.
1302 @note This operation destroys the current update operator (if owned). */
1304
1305 /// Free the GridFunction update operator (if any), to save memory.
1306 virtual void UpdatesFinished() { Th.Clear(); }
1307
1308 /** Return update counter, similar to Mesh::GetSequence(). Used by
1309 GridFunction to check if it is up to date with the space. */
1310 long GetSequence() const { return sequence; }
1311
1312 /// Return whether or not the space is discontinuous (L2)
1313 bool IsDGSpace() const
1314 {
1315 return dynamic_cast<const L2_FECollection*>(fec) != NULL;
1316 }
1317
1318 /** In variable order spaces on nonconforming (NC) meshes, this function
1319 controls whether strict conformity is enforced in cases where coarse
1320 edges/faces have higher polynomial order than their fine NC neighbors.
1321 In the default (strict) case, the coarse side polynomial order is
1322 reduced to that of the lowest order fine edge/face, so all fine
1323 neighbors can interpolate the coarse side exactly. If relaxed == true,
1324 some discontinuities in the solution in such cases are allowed and the
1325 coarse side is not restricted. For an example, see
1326 https://github.com/mfem/mfem/pull/1423#issuecomment-621340392 */
1327 void SetRelaxedHpConformity(bool relaxed = true)
1328 {
1329 relaxed_hp = relaxed;
1330 orders_changed = true; // force update
1331 Update(false);
1332 }
1333
1334 /// Save finite element space to output stream @a out.
1335 void Save(std::ostream &out) const;
1336
1337 /** @brief Read a FiniteElementSpace from a stream. The returned
1338 FiniteElementCollection is owned by the caller. */
1339 FiniteElementCollection *Load(Mesh *m, std::istream &input);
1340
1341 virtual ~FiniteElementSpace();
1342};
1343
1344/// @brief Return true if the mesh contains only one topology and the elements are tensor elements.
1346{
1347 Mesh & mesh = *fes.GetMesh();
1348 const bool mixed = mesh.GetNumGeometries(mesh.Dimension()) > 1;
1349 // Potential issue: empty local mesh --> no element 0.
1350 return !mixed &&
1351 dynamic_cast<const mfem::TensorBasisElement *>(fes.GetFE(0))!=nullptr;
1352}
1353
1354/// @brief Return LEXICOGRAPHIC if mesh contains only one topology and the elements are tensor
1355/// elements, otherwise, return NATIVE.
1356ElementDofOrdering GetEVectorOrdering(const FiniteElementSpace& fes);
1357
1358}
1359
1360#endif
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
Abstract base class BilinearFormIntegrator.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Rank 3 tensor (array of matrices)
Abstract base class that defines an interface for element restrictions.
A class that performs interpolation from a face E-vector to quadrature point values and/or derivative...
Base class for operators that extracts Face degrees of freedom.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition fe_coll.hpp:225
Derefinement operator, used by the friend class InterpolationGridTransfer.
Definition fespace.hpp:443
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition fespace.cpp:2037
DerefinementOperator(const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes, BilinearFormIntegrator *mass_integ)
TODO: Implement DofTransformation support.
Definition fespace.cpp:1936
GridFunction interpolation operator applicable after mesh refinement.
Definition fespace.hpp:417
virtual 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 fespace.cpp:1729
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, Table *old_elem_fos, int old_ndofs)
Definition fespace.cpp:1570
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition fespace.cpp:1665
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
void Save(std::ostream &out) const
Save finite element space to output stream out.
Definition fespace.cpp:3549
int GetEntityVDofs(int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID, int variant=0) const
Helper to get vertex, edge or face VDOFs (entity=0,1,2 resp.).
Definition fespace.cpp:976
void GetVDofs(int vd, Array< int > &dofs, int ndofs=-1) const
Returns the indices of all of the VDofs for the specified dimension 'vd'.
Definition fespace.cpp:193
DofTransformation DoFTrans
Definition fespace.hpp:275
int GetNVDofs() const
Number of all scalar vertex dofs.
Definition fespace.hpp:730
static int EncodeDof(int entity_base, int idx)
Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
Definition fespace.hpp:1013
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.cpp:1309
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition fespace.cpp:461
Array< char > var_face_orders
Definition fespace.hpp:259
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:581
void BuildNURBSFaceToDofTable() const
Generates partial face_dof table for a NURBS space.
Definition fespace.cpp:2325
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I, int skipfirst=0)
Definition fespace.cpp:807
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element,...
Definition fespace.hpp:1136
Array< int > dof_ldof_array
Definition fespace.hpp:268
int GetElementType(int i) const
Returns the type of element i.
Definition fespace.hpp:761
Array< StatelessDofTransformation * > DoFTransArray
Definition fespace.hpp:274
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Compute the full set of vdofs corresponding to each entry in dofs.
Definition fespace.cpp:213
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified edge.
Definition fespace.cpp:3149
Array< FaceQuadratureInterpolator * > E2BFQ_array
Definition fespace.hpp:311
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
Definition fespace.cpp:3204
Array< FaceQuadratureInterpolator * > E2IFQ_array
Definition fespace.hpp:310
NURBSExtension * GetNURBSext()
Definition fespace.hpp:562
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
Definition fespace.cpp:2846
int GetNumElementInteriorDofs(int i) const
Returns the number of degrees of freedom associated with the interior of the specified element.
Definition fespace.cpp:3119
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition fespace.hpp:769
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:773
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vert...
Definition fespace.cpp:322
NURBSExtension * NURBSext
Definition fespace.hpp:270
void GetUpdateOperator(OperatorHandle &T)
Return the update operator in the given OperatorHandle, T.
Definition fespace.hpp:1288
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.hpp:620
virtual int GetFaceDofs(int face, Array< int > &dofs, int variant=0) const
Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edge...
Definition fespace.cpp:2958
virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes,...
Definition fespace.cpp:3376
static void AdjustVDofs(Array< int > &vdofs)
Remove the orientation information encoded into an array of dofs Some basis function types have a rel...
Definition fespace.cpp:265
const Table & GetFaceToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each face in the m...
Definition fespace.hpp:1149
int GetEdgeOrder(int edge, int variant=0) const
Definition fespace.cpp:2699
bool Nonconforming() const
Definition fespace.hpp:566
void GetVertexVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified vertices.
Definition fespace.cpp:328
void GetVertexDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the specified vertices.
Definition fespace.cpp:3094
OperatorHandle L2E_lex
Definition fespace.hpp:293
void BuildFaceToDofTable() const
Definition fespace.cpp:426
int GetFaceOrder(int face, int variant=0) const
Returns the polynomial degree of the i'th face finite element.
Definition fespace.cpp:2710
int GetNEDofs() const
Number of all scalar edge-interior dofs.
Definition fespace.hpp:732
int GetNumBorderDofs(Geometry::Type geom, int order) const
Definition fespace.cpp:942
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Definition fespace.cpp:59
Array< int > dof_elem_array
Definition fespace.hpp:268
int FindFaceDof(int face, int ndof) const
Similar to FindEdgeDof, but used for mixed meshes too.
Definition fespace.hpp:373
static int MinOrder(VarOrderBits bits)
Return the minimum order (least significant bit set) in the bit mask.
Definition fespace.cpp:2503
int GetAttribute(int i) const
Definition fespace.hpp:785
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked...
Definition fespace.cpp:667
MFEM_DEPRECATED void RebuildElementToDofTable()
(
Definition fespace.cpp:452
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition fespace.hpp:710
bool orders_changed
True if at least one element order changed (variable-order space only).
Definition fespace.hpp:322
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos)
Calculate GridFunction restriction matrix after mesh derefinement.
Definition fespace.cpp:2095
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Returns the transformation defining the i-th element in the user-defined variable ElTr.
Definition fespace.hpp:778
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
Definition fespace.cpp:1528
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the interior of the specified edge.
Definition fespace.cpp:340
void GetTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes,...
Definition fespace.cpp:3349
SparseMatrix * D2C_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the continuous FE space of th...
Definition fespace.cpp:697
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition fespace.cpp:588
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition fespace.hpp:749
virtual void UpdateMeshPointer(Mesh *new_mesh)
Definition fespace.cpp:3544
const NURBSExtension * GetNURBSext() const
Definition fespace.hpp:561
int GetDegenerateFaceDofs(int index, Array< int > &dofs, Geometry::Type master_geom, int variant) const
Definition fespace.cpp:900
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition fespace.hpp:597
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition fespace.cpp:1404
int GetBdrAttribute(int i) const
Definition fespace.hpp:787
int GetEntityDofs(int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID, int variant=0) const
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition fespace.cpp:951
virtual const Operator * GetRestrictionOperator() const
An abstract operator that performs the same action as GetRestrictionMatrix.
Definition fespace.hpp:616
int GetLocalDofForDof(int i) const
Return the local dof index in the first element that contains dof i.
Definition fespace.hpp:1163
static constexpr int MaxVarOrder
Definition fespace.hpp:346
int GetElementForDof(int i) const
Return the index of the first element that contains dof i.
Definition fespace.hpp:1159
virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes, const Array< int > *perm)
Copies the prolongation and restriction matrices from fes.
Definition fespace.cpp:98
Array< char > var_edge_orders
Definition fespace.hpp:259
std::unique_ptr< Operator > R_transpose
Operator computing the action of the transpose of the restriction.
Definition fespace.hpp:287
bool Conforming() const
Definition fespace.hpp:565
void UpdateElementOrders()
Resize the elem_order array on mesh change.
Definition fespace.cpp:3417
void MakeVDimMatrix(SparseMatrix &mat) const
Replicate 'mat' in the vector dimension, according to vdim ordering mode.
Definition fespace.cpp:1274
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition fespace.hpp:746
SparseMatrix * H2L_GlobalRestrictionMatrix(FiniteElementSpace *lfes)
Construct the restriction matrix from the FE space given by (*this) to the lower degree FE space give...
Definition fespace.cpp:760
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE spa...
Definition fespace.cpp:729
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition fespace.hpp:231
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller.
Definition fespace.cpp:3618
int VDofToDof(int vdof) const
Compute the inverse of the Dof to VDof mapping for a single index vdof.
Definition fespace.hpp:995
void BuildBdrElementToDofTable() const
Definition fespace.cpp:387
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
Array< QuadratureInterpolator * > E2Q_array
Definition fespace.hpp:309
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:725
NURBSExtension * StealNURBSext()
Definition fespace.cpp:2290
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition fespace.hpp:782
void BuildDofToArrays()
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof...
Definition fespace.cpp:481
void GetPatchDofs(int patch, Array< int > &dofs) const
Returns indices of degrees of freedom for NURBS patch index patch. Cartesian ordering is used,...
Definition fespace.cpp:3161
FiniteElementSpace & operator=(const FiniteElementSpace &)=delete
Copy assignment not supported.
const Operator * GetRestrictionTransposeOperator() const
Return an operator that performs the transpose of GetRestrictionOperator.
Definition fespace.cpp:1324
Array< char > elem_order
Definition fespace.hpp:246
int FirstFaceDof(int face, int variant=0) const
Definition fespace.hpp:376
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition fespace.cpp:886
Array< int > face_to_be
Definition fespace.hpp:272
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
Definition fespace.cpp:2071
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition fespace.cpp:985
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition fespace.hpp:234
Table var_face_dofs
NOTE: also used for spaces with mixed faces.
Definition fespace.hpp:255
int GetOrder(int i) const
Returns the polynomial degree of the i'th finite element.
Definition fespace.hpp:696
OperatorHandle L2E_nat
The element restriction operators, see GetElementRestriction().
Definition fespace.hpp:293
FiniteElementSpace(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Construct a NURBS FE space based on the given NURBSExtension, ext.
Definition fespace.hpp:550
int GetConformingVSize() const
Definition fespace.hpp:722
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition fespace.hpp:1306
int * bdofs
internal DOFs of elements if mixed/var-order; NULL otherwise
Definition fespace.hpp:250
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition fespace.cpp:1336
const FiniteElement * GetEdgeElement(int i, int variant=0) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th edge in the ...
Definition fespace.cpp:3267
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
Definition fespace.hpp:1295
int GetEdgeDofs(int edge, Array< int > &dofs, int variant=0) const
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vert...
Definition fespace.cpp:3046
std::tuple< bool, ElementDofOrdering, FaceType, L2FaceValues > key_face
The face restriction operators, see GetFaceRestriction().
Definition fespace.hpp:295
std::unique_ptr< SparseMatrix > cR
Conforming restriction matrix such that cR.cP=I.
Definition fespace.hpp:282
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the interior of the specified element.
Definition fespace.cpp:334
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edge...
Definition fespace.cpp:316
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition fespace.cpp:3439
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition fespace.hpp:290
void CalcEdgeFaceVarOrders(Array< VarOrderBits > &edge_orders, Array< VarOrderBits > &face_orders) const
Definition fespace.cpp:2513
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos)
Definition fespace.cpp:1551
std::uint64_t VarOrderBits
Bit-mask representing a set of orders needed by an edge/face.
Definition fespace.hpp:345
void GetElementVertices(int i, Array< int > &vertices) const
Returns the vertices of element i.
Definition fespace.hpp:765
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const Table *coarse_elem_fos, const DenseTensor localP[]) const
Definition fespace.cpp:1464
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition fespace.hpp:242
std::unique_ptr< SparseMatrix > cP
Definition fespace.hpp:280
void SetUpdateOperatorType(Operator::Type tid)
Specify the Operator::Type to be used by the update operators.
Definition fespace.hpp:1303
void AddEdgeFaceDependencies(SparseMatrix &deps, Array< int > &master_dofs, const FiniteElement *master_fe, Array< int > &slave_dofs, int slave_face, const DenseMatrix *pm) const
Definition fespace.cpp:832
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
Definition fespace.cpp:3276
void SetRelaxedHpConformity(bool relaxed=true)
Definition fespace.hpp:1327
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Definition fespace.hpp:541
int MakeDofTable(int ent_dim, const Array< int > &entity_orders, Table &entity_dofs, Array< char > *var_ent_order)
Definition fespace.cpp:2629
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
Definition fespace.cpp:648
int GetNFDofs() const
Number of all scalar face-interior dofs.
Definition fespace.hpp:734
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
Definition fespace.cpp:176
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition fespace.hpp:228
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the ...
Definition fespace.cpp:3237
Ordering::Type ordering
Definition fespace.hpp:239
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
For a partially conforming FE space, convert a marker array (nonzero entries are true) on the partial...
Definition fespace.cpp:680
void SetElementOrder(int i, int p)
Sets the order of the i'th finite element.
Definition fespace.cpp:149
const SparseMatrix * GetHpConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.cpp:1317
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition fespace.hpp:757
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.cpp:1302
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:727
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetNV() const
Returns number of vertices in the mesh.
Definition fespace.hpp:737
int GetElementOrderImpl(int i) const
Return element order: internal version of GetElementOrder without checks.
Definition fespace.cpp:187
int GetNVariants(int entity, int index) const
Return number of possible DOF variants for edge/face (var. order spaces).
Definition fespace.cpp:2726
int GetMaxElementOrder() const
Return the maximum polynomial order.
Definition fespace.hpp:577
const Table & GetBdrElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each boundary mesh...
Definition fespace.hpp:1141
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
const Table * GetElementToFaceOrientationTable() const
Definition fespace.hpp:1132
DofTransformation * GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element 'bel'. The returned indices are offsets in...
Definition fespace.cpp:2950
void GetBoundaryTrueDofs(Array< int > &boundary_dofs, int component=-1)
Get a list of all boundary true dofs, boundary_dofs. For spaces with 'vdim' > 1, the 'component' para...
Definition fespace.cpp:632
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
const FaceQuadratureInterpolator * GetFaceQuadratureInterpolator(const IntegrationRule &ir, FaceType type) const
Return a FaceQuadratureInterpolator that interpolates E-vectors to quadrature point values and/or der...
Definition fespace.cpp:1433
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition fespace.hpp:1313
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified face.
Definition fespace.cpp:3125
virtual const SparseMatrix * GetHpRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition fespace.hpp:624
int GetNConformingDofs() const
Definition fespace.cpp:1330
std::unordered_map< const key_face, FaceRestriction *, key_hash > map_L2F
Definition fespace.hpp:306
static int DecodeDof(int dof)
Helper to return the DOF associated with a sign encoded DOF.
Definition fespace.hpp:1017
static int DecodeDof(int dof, real_t &sign)
Helper to determine the DOF and sign of a sign encoded DOF.
Definition fespace.hpp:1021
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified element.
Definition fespace.cpp:3104
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
Definition fespace.cpp:2198
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets int...
Definition fespace.cpp:303
void ConvertFromConformingVDofs(const Array< int > &cdofs, Array< int > &dofs)
For a partially conforming FE space, convert a marker array (nonzero entries are true) on the conform...
Definition fespace.cpp:688
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
Mark degrees of freedom associated with boundary elements with the specified boundary attributes (mar...
Definition fespace.cpp:514
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
Definition fespace.cpp:249
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
Return an Operator that converts L-vectors to E-vectors on each face.
Definition fespace.cpp:1369
const Operator * GetUpdateOperator()
Get the GridFunction update operator.
Definition fespace.hpp:1285
int FindEdgeDof(int edge, int ndof) const
Definition fespace.hpp:369
void GetPatchVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom in vdofs for NURBS patch i.
Definition fespace.cpp:310
void BuildElementToDofTable() const
Definition fespace.cpp:346
std::unique_ptr< SparseMatrix > cR_hp
A version of the conforming restriction matrix for variable-order spaces.
Definition fespace.hpp:284
int FindDofs(const Table &var_dof_table, int row, int ndof) const
Search row of a DOF table for a DOF set of size 'ndof', return first DOF.
Definition fespace.cpp:2682
Abstract class for all finite elements.
Definition fe_base.hpp:239
static const int NumGeom
Definition geom.hpp:42
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition transfer.hpp:124
A standard isoparametric element transformation.
Definition eltrans.hpp:363
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Definition lor.hpp:23
Mesh data type.
Definition mesh.hpp:56
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition mesh.cpp:7316
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition mesh.cpp:7321
bool Conforming() const
Definition mesh.hpp:2228
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition mesh.cpp:6250
int GetAttribute(int i) const
Return the attribute of element i.
Definition mesh.hpp:1333
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition mesh.hpp:1438
bool Nonconforming() const
Definition mesh.hpp:2229
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition mesh.hpp:1339
virtual int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type, does not count master nonconforming face...
Definition mesh.cpp:6266
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:357
ElementTransformation * GetBdrElementTransformation(int i)
Returns a pointer to the transformation defining the i-th boundary element.
Definition mesh.cpp:506
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition mesh.hpp:1223
void Swap(Mesh &other, bool non_geometry)
Definition mesh.cpp:10378
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1229
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition mesh.cpp:6961
Pointer to an Operator of a specified type.
Definition handle.hpp:34
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition handle.hpp:120
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition handle.hpp:132
Operator * Ptr() const
Access the underlying Operator pointer.
Definition handle.hpp:87
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged.
Definition handle.hpp:124
Abstract operator.
Definition operator.hpp:25
Type
Enumeration defining IDs for some classes derived from Operator.
Definition operator.hpp:284
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition fespace.hpp:30
static int Map(int ndofs, int vdim, int dof, int vd)
static void DofsToVDofs(int ndofs, int vdim, Array< int > &dofs)
Type
Ordering methods:
Definition fespace.hpp:34
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition transfer.hpp:429
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
Class representing the storage layout of a QuadratureFunction.
Definition qspace.hpp:120
Data type sparse matrix.
Definition sparsemat.hpp:51
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:187
Vector data type.
Definition vector.hpp:80
int index(int i, int j, int nx, int ny)
Definition life.cpp:235
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition fespace.hpp:1345
QVectorLayout
Type describing possible layouts for Q-vectors.
Definition fespace.hpp:53
@ byNODES
NQPT x VDIM x NE (values) / NQPT x VDIM x DIM x NE (grads)
@ byVDIM
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
float real_t
Definition config.hpp:43
ElementDofOrdering GetEVectorOrdering(const FiniteElementSpace &fes)
Return LEXICOGRAPHIC if mesh contains only one topology and the elements are tensor elements,...
Definition fespace.cpp:3719
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:75
@ LEXICOGRAPHIC
Lexicographic ordering for tensor-product FiniteElements.
@ NATIVE
Native ordering as defined by the FiniteElement.
FaceType
Definition mesh.hpp:47
real_t p(const Vector &x, real_t t)
std::size_t operator()(const key_face &k) const
Definition fespace.hpp:298