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