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