PTEMS  0.1.0-dev+git.81fd0e4
PolyTopic Element Method Solver
ptems::DiscreteFunctionSpaceInterface< DIM, X > Class Template Referenceabstract

An interface declaration for DiscreteFunctionSpace which removes the size of the Codomain. More...

#include <ptems/space.hpp>

Inheritance diagram for ptems::DiscreteFunctionSpaceInterface< DIM, X >:

Public Member Functions

virtual std::size_t NumberGlobalDoFs () const =0
 Gets the global number of degrees of freedom for the space. More...
 
virtual bool IsContinuous () const =0
 Gets if the functions space is continuous over element boundaries. More...
 
virtual bool CanIncrementPolynomialDegree () const
 Gets whether the space can be modified by incrementing/decrementing the polynomial degree uniformly for any component. More...
 
virtual bool CanIncrementPolynomialDegree ([[maybe_unused]] std::size_t component) const
 Gets whether the space can be modified by incrementing/decrementing the polynomial degree uniformly for the specified codomain component. More...
 
virtual bool SupportsVariablePolynomialDegree () const
 Gets whether the space can be modified with varying polynomial degree for each element for any component. More...
 
virtual bool SupportsVariablePolynomialDegree ([[maybe_unused]] std::size_t component) const
 Gets whether the space can be modified with varying polynomial degree for each element for the specified codomain component. More...
 

Protected Member Functions

virtual void ComputeStronglyImposedDoFs ([[maybe_unused]] std::unordered_map< std::size_t, X > &dofs, [[maybe_unused]] std::size_t offset) const
 Computes the set of global DoFs with strongly imposed Dirichlet boundary conditions. More...
 
virtual void ComputeGlobalDoFOffset (std::size_t *dofOffset, std::size_t offset) const =0
 Computes the offset onto the global degrees of freedom for each component. More...
 
virtual std::size_tComputeDoFMapping (std::size_t *dofMapping, std::size_t elementIdx, std::size_t offset) const =0
 Computes the mapping between local dof numbers for an element and the global dof numbers. More...
 
virtual X * ComputeBasisFunctions (X *basis, std::size_t elementIdx, const Vector< DIM > &pt)=0
 Computes the basis functions for each component at the specified point of the specified element. More...
 
virtual FuncAndGradData< DIM, X > * ComputeGradBasisFunctions (FuncAndGradData< DIM, X > *basis, std::size_t elementIdx, const Vector< DIM > &pt)=0
 Computes the basis functions and first order derivatives for each component at the specified point of the specified element. More...
 
virtual FuncAndGradData< DIM, X > * ComputeGradValueBasisFunctions (FuncAndGradData< DIM, X > *basis, std::size_t elementIdx, const Vector< DIM > &pt)
 Computes the basis functions and first order derivatives of those basis functions for each component at the specified point of the specified element. More...
 
virtual FuncGradAndHessianData< DIM, X > * ComputeHessianBasisFunctions ([[maybe_unused]] FuncGradAndHessianData< DIM, X > *basis, [[maybe_unused]] std::size_t elementIdx, [[maybe_unused]] const Vector< DIM > &pt)
 Computes the basis functions, first order derivatives, and second order derivatives for each component at the specified point of the specified element. More...
 
virtual std::vector< X > * ComputeDoFAtBasis (std::vector< X > *dofs, [[maybe_unused]] std::size_t elementIdx)
 Computes the basis functions of the space evaluated by the functional for each degree of freedom. More...
 
virtual std::vector< X > * ComputeDoFAtBasisProjection (std::vector< X > *dofs, std::size_t elementIdx)
 Computes the value projection of the basis functions of the space evaluated by the functional for each degree of freedom. More...
 
virtual Vector< DIM > * ComputeLocalAnalyticityEstimate (Vector< DIM > *analyticity, [[maybe_unused]] std::size_t elementIdx, [[maybe_unused]] const X *dofs)
 Computes the local analyticity estimate of the element using Legendre coefficients. More...
 
virtual std::size_t CodomainDimensionSize () const =0
 Gets ths size of the codomain of the fuctions in the space. More...
 
virtual std::size_tFillPolynomialDegree (std::size_t element, std::size_t *nextPolydeg) const =0
 Populates an array with polynomial degree for each component of the result for the specified element. More...
 
virtual void SetPolynomialDegree ([[maybe_unused]] const std::map< std::size_t, const std::size_t * > &elementPolydegs, const DoFChangeList< X > &dofs)
 Set the polynomial degree of the specified elements and update the DoFs of the specified functions. More...
 
virtual void SetPolynomialDegree ([[maybe_unused]] const std::size_t *polydeg, const DoFChangeList< X > &dofs)
 Set the polynomial degree of the all elements and update the DoFs of the specified functions. More...
 
virtual void ComputeNumberLocalDoFs (std::size_t elementIdx, std::size_t *numberDofs) const =0
 Gets the number of degrees of freedom for an element in the space per component. More...
 
virtual std::vector< X >::iterator Interpolate (typename std::vector< X >::iterator dofsBegin, typename std::vector< X >::iterator dofsEnd, const X *value) const =0
 Fill the degrees of freedom for this space by interpolating the specified constant function. More...
 
virtual void OnMeshChanged (const PFEMesh< DIM > &previousMesh, const typename FEMesh< DIM >::Modifications &changes, const DoFChangeList< X > &dofs)=0
 Event notification from mesh that the mesh has been changed. More...
 

Friends

template<std::size_t D, typename T , std::size_t N>
class DiscreteCartesianProductSpace
 

Detailed Description

template<std::size_t DIM, typename X>
class ptems::DiscreteFunctionSpaceInterface< DIM, X >

An interface declaration for DiscreteFunctionSpace which removes the size of the Codomain.

Template Parameters
DIMThe dimension of the domain (and mesh)
XThe type of the result of the functions (double, std::complex, etc)

Member Function Documentation

◆ CanIncrementPolynomialDegree() [1/2]

template<std::size_t DIM, typename X >
virtual bool ptems::DiscreteFunctionSpaceInterface< DIM, X >::CanIncrementPolynomialDegree ( ) const
inlinevirtual

Gets whether the space can be modified by incrementing/decrementing the polynomial degree uniformly for any component.

Note
Must return true if SupportsVariablePolynomialDegree() returns true
Returns
true if the polynomial degree can be uniformly incremented; false otherwise

Reimplemented in ptems::DiscreteCartesianProductSpace< DIM, X, N >, ptems::VariableDegreeDiscreteFunctionSpace< DIM, X, N >, ptems::VariableDegreeDiscreteFunctionSpace< 1, X, N >, ptems::VariableDegreeDiscreteFunctionSpace< 2, X, N >, and ptems::VariableDegreeDiscreteFunctionSpace< DIM, double, 1 >.

◆ CanIncrementPolynomialDegree() [2/2]

template<std::size_t DIM, typename X >
virtual bool ptems::DiscreteFunctionSpaceInterface< DIM, X >::CanIncrementPolynomialDegree ( [[maybe_unused] ] std::size_t  component) const
inlinevirtual

Gets whether the space can be modified by incrementing/decrementing the polynomial degree uniformly for the specified codomain component.

Parameters
componentIndex of component to check
Returns
true if the polynomial degree of the specified component can be uniformly incremented; false otherwise

◆ CodomainDimensionSize()

template<std::size_t DIM, typename X >
virtual std::size_t ptems::DiscreteFunctionSpaceInterface< DIM, X >::CodomainDimensionSize ( ) const
protectedpure virtual

Gets ths size of the codomain of the fuctions in the space.

Implemented in ptems::DiscreteFunctionSpace< DIM, X, N >, and ptems::DiscreteFunctionSpace< DIM, double, 1 >.

◆ ComputeBasisFunctions()

template<std::size_t DIM, typename X >
virtual X* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeBasisFunctions ( X *  basis,
std::size_t  elementIdx,
const Vector< DIM > &  pt 
)
protectedpure virtual

Computes the basis functions for each component at the specified point of the specified element.

The point is the local point defined on the reference element (for elements with a single reference element) or the bounding box (for elements with multiple reference elements).

The passed pointer is a continuous array of all the basis functions for the first component, followed by all the basis functions for the second etc. The number of basis functions MUST be the same as ComputeNumberLocalDoFs calculates

Warning
Subclasses MUST add the number of basis given by ComputeNumberLocalDoFs for each component and return basis plus the sum of the number of basis functions of all components as computed by ComputeNumberLocalDoFs. Failure to do this is undefined behaviour (and could lead to segfaults)
Note
For virtual element the returned values are the value projections of the basis functions
Parameters
basisPointer to the first basis function to compute. This pointer can be treated as if it is an array containing the number of DoFs for the first component as computed by ComputeNumberLocalDoFs(), followed by the number of DoFs for the second component as computed by ComputeNumberLocalDoFs(), etc.
elementIdxThe element index to get the basis functions for
ptThe local point in the reference element or bounding box to get the basis functions for
Returns
Pointer to past the end of the array of basis functions computed. Must be equal to basis+x, where x is the sum of the array entries computed by ComputeNumberLocalDoFs

Implemented in ptems::DiscontinuousSpace< Polynomials, DIM, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ ComputeDoFAtBasis()

template<std::size_t DIM, typename X >
virtual std::vector<X>* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeDoFAtBasis ( std::vector< X > *  dofs,
[[maybe_unused] ] std::size_t  elementIdx 
)
inlineprotectedvirtual

Computes the basis functions of the space evaluated by the functional for each degree of freedom.

Parameters
dofsPointer to the vector to store each basis for the first DoF functional to compute. This pointer can be treated as if it is an array containing the number of DoFs for the first component as computed by ComputeNumberLocalDoFs(), followed by the number of DoFs for the second component as computed by ComputeNumberLocalDoFs(), etc.. Each vector will already be sized to same size as ComputeNumberLocalDoFs() for the component and filled wth zero
elementIdxThe element index to get the basis functions/dofs for
Returns
Pointer to past the end of the array of basis functions computed. Must be equal to dofs+x, where x is the sum of the array entries computed by ComputeNumberLocalDoFs

◆ ComputeDoFAtBasisProjection()

template<std::size_t DIM, typename X >
virtual std::vector<X>* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeDoFAtBasisProjection ( std::vector< X > *  dofs,
std::size_t  elementIdx 
)
inlineprotectedvirtual

Computes the value projection of the basis functions of the space evaluated by the functional for each degree of freedom.

Note
This function is only relevant for virtual element spaces. For finite element spaces it simply returns EvaluateDoFAtBasis();
Parameters
dofsPointer to the vector to store each basis for the first DoF functional to compute. This pointer can be treated as if it is an array containing the number of DoFs for the first component as computed by ComputeNumberLocalDoFs(), followed by the number of DoFs for the second component as computed by ComputeNumberLocalDoFs(), etc. Each vector will already be sized to same size as ComputeNumberLocalDoFs() for the component and filled wth zero
elementIdxThe element index to get the basis functions/dofs for
Returns
Pointer to past the end of the array of basis functions computed. Must be equal to dofs+x, where x is the sum of the array entries computed by ComputeNumberLocalDoFs

Reimplemented in ptems::VirtualElementSpace< HBasis, 2, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ ComputeDoFMapping()

template<std::size_t DIM, typename X >
virtual std::size_t* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeDoFMapping ( std::size_t dofMapping,
std::size_t  elementIdx,
std::size_t  offset 
) const
protectedpure virtual

Computes the mapping between local dof numbers for an element and the global dof numbers.

The passed pointer is a continuous array of all the dofs for the first component, followed by all the dofs for the second etc. The number of dofs MUST be the same as ComputeNumberLocalDoFs calculates

Note
Implementations of this function can cache this value for performance IF the number only changes when the mesh changes. The implementation must then invalidate the cached value when OnMeshChanged is called.
Exceptions
std::logic_errorIf the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space.
std::out_of_rangeif elementIdx is not less than the number of elements in the underlying mesh
Parameters
dofMappingPointer to the first global dof index to compute. This pointer can be treated as if it is an array containing the number of DoFs for the first component as computed by ComputeNumberLocalDoFs(), followed by the number of DoFs for the second component as computed by ComputeNumberLocalDoFs(), etc.
elementIdxIndex of the element to get the local to global mapping of the degrees of freedom for
offsetOffset to add to the global DoF numbers
Returns
Pointer to past the end of the array of global dof indices computed. Must be equal to dofMapping+x, where x is the sum of the array entries computed by ComputeNumberLocalDoFs

Implemented in ptems::VirtualElementSpace< HBasis, 2, X, N >, ptems::PiecewiseConstant< DIM, X, N >, ptems::DiscontinuousSpace< Polynomials, DIM, X, N >, ptems::DiscontinuousSpace< LegendrePolynomials< DIM >, DIM, double, 1 >, ptems::DiscontinuousSpace< LagrangePolynomials< DIM >, DIM, double, 1 >, ptems::ContinuousLagrangeSpace< 2, X, N >, ptems::ContinuousLagrangeSpace< 1, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ ComputeGlobalDoFOffset()

template<std::size_t DIM, typename X >
virtual void ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeGlobalDoFOffset ( std::size_t dofOffset,
std::size_t  offset 
) const
protectedpure virtual

Computes the offset onto the global degrees of freedom for each component.

Note
Implementations of this function can cache this value for performance IF the number only changes when the mesh changes. The implementation must then invalidate the cached value when OnMeshChanged is called.
Exceptions
std::logic_errorIf the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space.
Parameters
dofOffsetPointer to array of length of number of components + 1 to store the the offset into the global DoFs for the DoFs for each component (first N entries) and the total number of globals DoFs (last entry)
offsetOffset to add to the global DoF numbers

Implemented in ptems::VirtualElementSpace< HBasis, 2, X, N >, ptems::PiecewiseConstant< DIM, X, N >, ptems::DiscontinuousSpace< Polynomials, DIM, X, N >, ptems::DiscontinuousSpace< LegendrePolynomials< DIM >, DIM, double, 1 >, ptems::DiscontinuousSpace< LagrangePolynomials< DIM >, DIM, double, 1 >, ptems::ContinuousLagrangeSpace< 2, X, N >, ptems::ContinuousLagrangeSpace< 1, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ ComputeGradBasisFunctions()

template<std::size_t DIM, typename X >
virtual FuncAndGradData<DIM,X>* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeGradBasisFunctions ( FuncAndGradData< DIM, X > *  basis,
std::size_t  elementIdx,
const Vector< DIM > &  pt 
)
protectedpure virtual

Computes the basis functions and first order derivatives for each component at the specified point of the specified element.

The point is the local point defined on the reference element (for elements with a single reference element) or the bounding box (for elements with multiple reference elements).

The passed pointer is a continuous array of all the basis functions for the first component, followed by all the basis functions for the second etc. The number of basis functions MUST be the same as ComputeNumberLocalDoFs calculates

Note
For virtual element the returned values are the value and gradient projections of the basis functions
Warning
Subclasses MUST add the number of basis given by ComputeNumberLocalDoFs for each component and return basis plus the sum of the number of basis functions of all components as computed by ComputeNumberLocalDoFs. Failure to do this is undefined behaviour (and could lead to segfaults)
Parameters
basisPointer to the first basis function to compute. This pointer can be treated as if it is an array containing the number of DoFs for the first component as computed by ComputeNumberLocalDoFs(), followed by the number of DoFs for the second component as computed by ComputeNumberLocalDoFs(), etc.
elementIdxThe element index to get the basis functions for
ptThe local point in the reference element or bounding box to get the basis functions for
Returns
Pointer to past the end of the array of basis functions computed. Must be equal to basis+x, where x is the sum of the array entries computed by ComputeNumberLocalDoFs

Implemented in ptems::DiscontinuousSpace< Polynomials, DIM, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ ComputeGradValueBasisFunctions()

template<std::size_t DIM, typename X >
virtual FuncAndGradData<DIM,X>* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeGradValueBasisFunctions ( FuncAndGradData< DIM, X > *  basis,
std::size_t  elementIdx,
const Vector< DIM > &  pt 
)
inlineprotectedvirtual

Computes the basis functions and first order derivatives of those basis functions for each component at the specified point of the specified element.

The point is the local point defined on the reference element (for elements with a single reference element) or the bounding box (for elements with multiple reference elements).

The passed pointer is a continuous array of all the basis functions for the first component, followed by all the basis functions for the second etc. The number of basis functions MUST be the same as ComputeNumberLocalDoFs calculates

Note
Only relevant for virtual element spaces - the returned values are the value and gradient of the value projection of the basis functions; for other methods this returns the same as ComputeGradBasisFunctions
Warning
Subclasses MUST add the number of basis given by ComputeNumberLocalDoFs for each component and return basis plus the sum of the number of basis functions of all components as computed by ComputeNumberLocalDoFs. Failure to do this is undefined behaviour (and could lead to segfaults)
Parameters
basisPointer to the first basis function to compute. This pointer can be treated as if it is an array containing the number of DoFs for the first component as computed by ComputeNumberLocalDoFs(), followed by the number of DoFs for the second component as computed by ComputeNumberLocalDoFs(), etc.
elementIdxThe element index to get the basis functions for
ptThe local point in the reference element or bounding box to get the basis functions for
Returns
Pointer to past the end of the array of basis functions computed. Must be equal to basis+x, where x is the sum of the array entries computed by ComputeNumberLocalDoFs

◆ ComputeHessianBasisFunctions()

template<std::size_t DIM, typename X >
virtual FuncGradAndHessianData<DIM,X>* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeHessianBasisFunctions ( [[maybe_unused] ] FuncGradAndHessianData< DIM, X > *  basis,
[[maybe_unused] ] std::size_t  elementIdx,
[[maybe_unused] ] const Vector< DIM > &  pt 
)
inlineprotectedvirtual

Computes the basis functions, first order derivatives, and second order derivatives for each component at the specified point of the specified element.

The point is the local point defined on the reference element (for elements with a single reference element), or the bounding box (for elements with multiple reference elements).

The passed pointer is a continuous array of all the basis functions for the first component, followed by all the basis functions for the second etc. The number of basis functions MUST be the same as ComputeNumberLocalDoFs calculates

Warning
Subclasses MUST add the number of basis given by ComputeNumberLocalDoFs for each component and return basis plus the sum of the number of basis functions of all components as computed by ComputeNumberLocalDoFs. Failure to do this is undefined behaviour (and could lead to segfaults)
Exceptions
std::logic_errorIf second order derivatives are not supported by the space
Parameters
basisPointer to the first basis function to compute. This pointer can be treated as if it is an array containing the number of DoFs for the first component as computed by ComputeNumberLocalDoFs(), followed by the number of DoFs for the second component as computed by ComputeNumberLocalDoFs(), etc.
elementIdxThe element index to get the basis functions for
ptThe local point in the reference element or bounding box to get the basis functions for
Returns
Pointer to past the end of the array of basis functions computed. Must be equal to basis+x, where x is the sum of the array entries computed by ComputeNumberLocalDoFs

◆ ComputeLocalAnalyticityEstimate()

template<std::size_t DIM, typename X >
virtual Vector<DIM>* ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeLocalAnalyticityEstimate ( Vector< DIM > *  analyticity,
[[maybe_unused] ] std::size_t  elementIdx,
[[maybe_unused] ] const X *  dofs 
)
inlineprotectedvirtual

Computes the local analyticity estimate of the element using Legendre coefficients.

Parameters
elementIdxElement to estimate analyticity
analyticityPointer to the array to store analyticity estimate, per variable, in range \((0,1]\) (or negative if analyticity cannot be computed)
dofsPointer to array containing degrees of freedom
Returns
Pointer to past the end of the array of analyticities computed

◆ ComputeNumberLocalDoFs()

template<std::size_t DIM, typename X >
virtual void ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeNumberLocalDoFs ( std::size_t  elementIdx,
std::size_t numberDofs 
) const
protectedpure virtual

Gets the number of degrees of freedom for an element in the space per component.

Note
Implementations of this function can cache this value for performance IF the number only changes when the mesh changes. The implementation must then invalidate the cached value when OnMeshChanged is called.
Exceptions
std::logic_errorIf the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space.
std::out_of_rangeif elementIdx is not less than the number of elements in the underlying mesh
Parameters
elementIdxIndex of the element to get the number of local degrees of freedom for
numberDofsPointer to array of length of number of components to store the local number of degrees of freedom for the specified element per component

Implemented in ptems::PiecewiseConstant< DIM, X, N >, ptems::VirtualElementSpace< HBasis, 2, X, N >, ptems::DiscontinuousSpace< Polynomials, DIM, X, N >, ptems::DiscontinuousSpace< LegendrePolynomials< DIM >, DIM, double, 1 >, ptems::DiscontinuousSpace< LagrangePolynomials< DIM >, DIM, double, 1 >, ptems::ContinuousLagrangeSpace< 2, X, N >, ptems::ContinuousLagrangeSpace< 1, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ ComputeStronglyImposedDoFs()

template<std::size_t DIM, typename X >
virtual void ptems::DiscreteFunctionSpaceInterface< DIM, X >::ComputeStronglyImposedDoFs ( [[maybe_unused] ] std::unordered_map< std::size_t, X > &  dofs,
[[maybe_unused] ] std::size_t  offset 
) const
inlineprotectedvirtual

Computes the set of global DoFs with strongly imposed Dirichlet boundary conditions.

Parameters
dofsSet to update with global dof number and value for strongly imposed Dirichlet BC
offsetOffset to add to the global DoF numbers

◆ FillPolynomialDegree()

template<std::size_t DIM, typename X >
virtual std::size_t* ptems::DiscreteFunctionSpaceInterface< DIM, X >::FillPolynomialDegree ( std::size_t  element,
std::size_t nextPolydeg 
) const
protectedpure virtual

Populates an array with polynomial degree for each component of the result for the specified element.

Parameters
elementThe index of the element to get the polynomial degree for
nextPolydegPoint into an array to store the polynomial degree. The array will have enough space to store CodomainDimensionSize() values
Returns
std::size_t* The pointer to the array entry AFTER the last polynomial degree entered. Must be equal to nextPolydeg+CodomainDimensionSize()

Implemented in ptems::DiscreteFunctionSpace< DIM, X, N >, and ptems::DiscreteFunctionSpace< DIM, double, 1 >.

◆ Interpolate()

template<std::size_t DIM, typename X >
virtual std::vector<X>::iterator ptems::DiscreteFunctionSpaceInterface< DIM, X >::Interpolate ( typename std::vector< X >::iterator  dofsBegin,
typename std::vector< X >::iterator  dofsEnd,
const X *  value 
) const
protectedpure virtual

Fill the degrees of freedom for this space by interpolating the specified constant function.

Exceptions
std::logic_errorIf dofsEnd-dofsBegin is less than the number of degrees of freedom to be populated
Warning
Note that dofsEnd-dofsBegin may be larger than the number of degrees of freedom in the space. Implementations of this function MUST ONLY populate the same number degrees of degrees of freedom as NumberGlobalDoFs() and return the iterator to dofsBegin+NumberGlobalDoFs()
Parameters
dofsBeginIterator pointing to the start of the degrees of freedom to fill
dofsEndIterator pointing to the end of the degrees of freedom to fill
valuePointer to array containing the constants for each component
Returns
Iterator pointing to past the last degree of freedom populated by the interpolation. Must be dofsBegin+NumberGlobalDoFs().

Implemented in ptems::VirtualElementSpace< HBasis, 2, X, N >, ptems::DiscontinuousSpace< Polynomials, DIM, X, N >, ptems::ContinuousLagrangeSpace< 2, X, N >, ptems::ContinuousLagrangeSpace< 1, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ IsContinuous()

◆ NumberGlobalDoFs()

template<std::size_t DIM, typename X >
virtual std::size_t ptems::DiscreteFunctionSpaceInterface< DIM, X >::NumberGlobalDoFs ( ) const
pure virtual

Gets the global number of degrees of freedom for the space.

Note
Implementations of this function can cache this value for performance IF the number only changes when the mesh changes. The implementation must then invalidate the cached value when OnMeshChanged is called.
Exceptions
std::logic_errorIf the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space.
Returns
Global number of degrees of freedom

Implemented in ptems::VirtualElementSpace< HBasis, 2, X, N >, ptems::PiecewiseConstant< DIM, X, N >, ptems::DiscontinuousSpace< Polynomials, DIM, X, N >, ptems::DiscontinuousSpace< LegendrePolynomials< DIM >, DIM, double, 1 >, ptems::DiscontinuousSpace< LagrangePolynomials< DIM >, DIM, double, 1 >, ptems::ContinuousLagrangeSpace< 2, X, N >, ptems::ContinuousLagrangeSpace< 1, X, N >, and ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ OnMeshChanged()

template<std::size_t DIM, typename X >
virtual void ptems::DiscreteFunctionSpaceInterface< DIM, X >::OnMeshChanged ( const PFEMesh< DIM > &  previousMesh,
const typename FEMesh< DIM >::Modifications &  changes,
const DoFChangeList< X > &  dofs 
)
protectedpure virtual

Event notification from mesh that the mesh has been changed.

The space should perform any necessary changes to conform to the new mesh.

The m_mesh variable will contain the new mesh.

Warning
The target vector in dofs may already contain values, this MUST NOT be overwritten, new values should instead be appended.
Exceptions
std::logic_errorIf the new mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space.
Parameters
previousMeshThe mesh before refinement
changesA FEMesh<DIM>::Modifications structure defining the change to the mesh
dofsMap (actually array of pairs) of the degree of freedom from a function which is part of this space before the change, to a (reference to a) vector to APPEND with the new degrees of freedom after mesh change

Implemented in ptems::DiscreteCartesianProductSpace< DIM, X, N >.

◆ SetPolynomialDegree() [1/2]

template<std::size_t DIM, typename X >
virtual void ptems::DiscreteFunctionSpaceInterface< DIM, X >::SetPolynomialDegree ( [[maybe_unused] ] const std::map< std::size_t, const std::size_t * > &  elementPolydegs,
const DoFChangeList< X > &  dofs 
)
inlineprotectedvirtual

Set the polynomial degree of the specified elements and update the DoFs of the specified functions.

Warning
The target vector in dofs may already contain values, this MUST NOT be overwritten, new values should instead be appended.
Parameters
elementPolydegs(Ordered) map of element number to array (of length N) of new polynomial degrees for each component
dofsMap (actually array of pairs) of the degree of freedom from a function which is part of this space before the change, to a (reference to a) vector to APPEND with the new degrees of freedom after mesh change

◆ SetPolynomialDegree() [2/2]

template<std::size_t DIM, typename X >
virtual void ptems::DiscreteFunctionSpaceInterface< DIM, X >::SetPolynomialDegree ( [[maybe_unused] ] const std::size_t polydeg,
const DoFChangeList< X > &  dofs 
)
inlineprotectedvirtual

Set the polynomial degree of the all elements and update the DoFs of the specified functions.

Warning
The target vector in dofs may already contain values, this MUST NOT be overwritten, new values should instead be appended.
Parameters
polydegArray (of length N) of new polynomial degrees for each component
dofsMap (actually array of pairs) of the degree of freedom from a function which is part of this space before the change, to a (reference to a) vector to APPEND with the new degrees of freedom after mesh change

◆ SupportsVariablePolynomialDegree() [1/2]

template<std::size_t DIM, typename X >
virtual bool ptems::DiscreteFunctionSpaceInterface< DIM, X >::SupportsVariablePolynomialDegree ( ) const
inlinevirtual

Gets whether the space can be modified with varying polynomial degree for each element for any component.

Note
If this returns true then CanIncrementPolynomialDegree must also return true
Returns
true if the polynomial degree can be different on different elements; false otherwise

Reimplemented in ptems::ContinuousLagrangeSpace< 2, X, N >, ptems::ContinuousLagrangeSpace< 1, X, N >, ptems::DiscreteCartesianProductSpace< DIM, X, N >, ptems::VariableDegreeDiscreteFunctionSpace< DIM, X, N >, ptems::VariableDegreeDiscreteFunctionSpace< 1, X, N >, ptems::VariableDegreeDiscreteFunctionSpace< 2, X, N >, and ptems::VariableDegreeDiscreteFunctionSpace< DIM, double, 1 >.

◆ SupportsVariablePolynomialDegree() [2/2]

template<std::size_t DIM, typename X >
virtual bool ptems::DiscreteFunctionSpaceInterface< DIM, X >::SupportsVariablePolynomialDegree ( [[maybe_unused] ] std::size_t  component) const
inlinevirtual

Gets whether the space can be modified with varying polynomial degree for each element for the specified codomain component.

Parameters
componentIndex of component to check
Returns
true if the polynomial degree of the specified component can be different on different elements; false otherwise

Reimplemented in ptems::DiscreteCartesianProductSpace< DIM, X, N >.


The documentation for this class was generated from the following file: