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

Defines a space of discrete functions defined as the Cartesian product of discrete spaces. More...

#include <ptems/space.hpp>

Inheritance diagram for ptems::DiscreteCartesianProductSpace< DIM, X, N >:
Collaboration diagram for ptems::DiscreteCartesianProductSpace< DIM, X, N >:

Public Types

typedef X CodomainType
 The type of the result of the functions. More...
 
typedef DiscreteFunction< DIM, X, N > FunctionType
 The type of the functions in this space. More...
 

Public Member Functions

virtual std::size_t NumberGlobalDoFs () const override
 Gets the global number of degrees of freedom for the space. More...
 
const std::array< std::size_t, N > & PolynomialDegree (std::size_t element) const override
 Gets the polynomial degree of the function space for the specified element. More...
 
bool IsContinuous () const override
 Gets if the functions space is continuous over element boundaries. More...
 
bool CanIncrementPolynomialDegree () const override
 Gets whether the space can be modified by incrementing/decrementing the polynomial degree uniformly for any component. More...
 
bool CanIncrementPolynomialDegree (std::size_t component) const override
 
bool SupportsVariablePolynomialDegree () const override
 Gets whether the space can be modified with varying polynomial degree for each element for any component. More...
 
bool SupportsVariablePolynomialDegree ([[maybe_unused]] std::size_t component) const override
 Gets whether the space can be modified with varying polynomial degree for each element for the specified codomain component. More...
 
const PFEMesh< DIM > & Mesh () const
 Gets the mesh the function space is over. More...
 
void AttachMeshListener ()
 Attach the listener for this space to the mesh. More...
 
virtual void PolynomialDegree (std::size_t element, const std::array< std::size_t, N > &polydeg)
 Sets the polynomial degree of the function space for the specified element. More...
 
void PolynomialDegree (std::size_t element, std::size_t polydeg)
 Sets the polynomial degree of the function space for the specified element. More...
 
template<typename... T, typename = std::enable_if_t<sizeof...(T) == N && (sizeof...(T) > 1)>>
void PolynomialDegree (std::size_t element, T... polydeg)
 Sets the polynomial degree of the function space for the specified element. More...
 
virtual void UniformPolynomialDegree (const std::array< std::size_t, N > &polydeg)
 Sets the polynomial degree of the function space for all elements. More...
 
void UniformPolynomialDegree (std::size_t polydeg)
 Sets the polynomial degree of the function space for all elements. More...
 
template<typename... T, typename = std::enable_if_t<sizeof...(T) == N && (sizeof...(T) > 1)>>
void UniformPolynomialDegree (T... polydeg)
 Sets the polynomial degree of the function space for all elements. More...
 
template<typename Func , typename = std::enable_if_t< std::is_same_v< decltype(std::declval<Func>()(std::declval<std::size_t>(), std::declval<std::array<std::size_t,N>&>())), bool>>
void Adapt (Func computePolydeg)
 Performs p-refinement on the space using a function which updates the polynomial degree of each element to the new polynomial degree. More...
 
template<typename Func , typename = std::enable_if_t<std::is_same_v<decltype(std::declval<Func>()(std::declval<std::size_t>())), AdaptationType>>
void Adapt (Func flagFunction, std::size_t minimum=1)
 Performs p-refinement on the space using a function for deciding which elements to refine. More...
 
template<typename Container , typename = typename Container::value_type>
void Adapt (const Container &flags, std::size_t minimum=1)
 Performs p-refinement on the space using a container of flags indicating whether to perform refinement (increment polynomial degree by one) or coarsen (decrement polynomial degree by one) More...
 
template<typename Container >
void Adapt (const Container &refine, const Container &coarsen, std::size_t minimum=1)
 Performs p-refinement on the space using a set of elements to refine (increment polynomial degree by one) and coarsen (decrement polynomial degree by one) More...
 
SpacePolynomialDegree PolynomialDegreeFunction (std::size_t component=0) const
 Gets the polynomial degree of the function space for the specified component. More...
 
DiscreteFunction< DIM, X, N > Function ()
 Creates a function in this function space. More...
 
template<typename... T, typename = std::enable_if_t<(sizeof...(T) == N) || (sizeof...(T) == 1)>>
DiscreteFunction< DIM, X, N > Interpolant (T... values)
 Creates a function in this function space, with initial constant function value. More...
 
void MeshChanged (const PFEMesh< DIM > &previousMesh, const typename FEMesh< DIM >::Modifications &changes) override
 Event notification from mesh that the mesh has been changed. More...
 
virtual void MeshChanged (const std::shared_ptr< FEMesh< DIM >> &previousMesh, const typename FEMesh< DIM >::Modifications &changes)=0
 Event notification from mesh that the mesh has been changed. More...
 
std::unordered_map< std::size_t, X > StronglyImposedDoFs () const
 Computes the set of global DoFs with strongly imposed Dirichlet boundary conditions. More...
 
std::array< std::size_t, N+1 > GlobalDoFOffset () const
 Gets the offset into the global DoFs for each component. More...
 
std::array< std::size_t, N > NumberLocalDoFs (std::size_t elementIdx) const
 Gets the number of degrees of freedom for an element in the space per component. More...
 
DoFData< N, std::size_tDoFMapping (std::size_t elementIdx) const
 Gets the mapping between local dof numbers for an element and the global dof numbers. More...
 
DoFData< N, X > BasisFunctions (std::size_t elementIdx, const Vector< DIM > &pt)
 Gets the values of the basis functions for each component at the specified point of the specified element. More...
 
DoFData< N, FuncAndGradData< DIM, X > > GradBasisFunctions (std::size_t elementIdx, const Vector< DIM > &pt)
 Gets the derivative of the basis functions in each coordinate direction for each component at the specified point in the specified element. More...
 
DoFData< N, FuncAndGradData< DIM, X > > GradValueBasisFunctions (std::size_t elementIdx, const Vector< DIM > &pt)
 Gets the derivative of the value of basis functions in each coordinate direction for each component at the specified point in the specified element. More...
 
DoFData< N, FuncGradAndHessianData< DIM, X > > HessianBasisFunctions (std::size_t elementIdx, const Vector< DIM > &pt)
 Gets the second order derivatives of the basis functions for each component at the specified point in the specified element. More...
 
DoFData< N, std::vector< X > > EvaluateDoFAtBasis (std::size_t elementIdx)
 Gets the basis functions of the space evaluated by the functional for each degree of freedom. More...
 
DoFData< N, std::vector< X > > EvaluateDoFAtBasisProjection (std::size_t elementIdx)
 Gets the value projection of the basis functions of the space evaluated by the functional for each degree of freedom. More...
 
template<std::size_t M>
std::array< Vector< DIM >, N > EstimateLocalAnalyticity (std::size_t elementIdx, const DiscreteFunction< DIM, X, M, N > &func)
 Estimates the local analyticity of the element using Legendre coefficients. 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...
 
shared_from_this (T... args)
 
weak_from_this (T... args)
 

Static Public Attributes

static constexpr std::size_t DomainDimension = DIM
 The dimension of the spacial domain the functions are over. More...
 
static constexpr std::size_t CodomainDimension = N
 The size of the resulting vector space for the functions. More...
 

Protected Types

typedef std::weak_ptr< detail::DiscreteFunctionImpl< DIM, X, N > > WeakFunction
 Weak pointer type to a function constructed from this space. More...
 
typedef std::set< WeakFunction, std::owner_less< WeakFunction > > FunctionList
 List type (actually set) of weak pointers to functions constructed from this space. More...
 

Protected Member Functions

virtual void ComputeStronglyImposedDoFs (std::unordered_map< std::size_t, X > &dofs, std::size_t offset) const override
 
virtual void ComputeGlobalDoFOffset (std::size_t *dofOffset, std::size_t offset) const override
 Computes the offset onto the global degrees of freedom for each component. More...
 
virtual void ComputeNumberLocalDoFs (std::size_t elementIdx, std::size_t *numberDofs) const override
 Gets the number of degrees of freedom for an element in the space per component. More...
 
virtual std::size_tComputeDoFMapping (std::size_t *dofMapping, std::size_t elementIdx, std::size_t offset) const override
 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) override
 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) override
 Computes the basis functions and first order derivatives for each component at the specified point of the specified element. More...
 
virtual FuncGradAndHessianData< DIM, X > * ComputeHessianBasisFunctions (FuncGradAndHessianData< DIM, X > *basis, std::size_t elementIdx, const Vector< DIM > &pt) override
 
virtual std::vector< X > * ComputeDoFAtBasis (std::vector< X > *dofs, std::size_t elementIdx) override
 
virtual std::vector< X > * ComputeDoFAtBasisProjection (std::vector< X > *dofs, std::size_t elementIdx) override
 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, std::size_t elementIdx, const X *dofs) override
 
virtual std::vector< X >::iterator Interpolate (typename std::vector< X >::iterator dofsBegin, typename std::vector< X >::iterator dofsEnd, const X *value) const override
 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 &elementMapping, const DoFChangeList< X > &dofs) override
 Event notification from mesh that the mesh has been changed. More...
 
virtual void SetPolynomialDegree (const std::map< std::size_t, const std::size_t * > &elementPolydegs, const DoFChangeList< X > &dofs) override
 
virtual void SetPolynomialDegree (const std::size_t *polydeg, const DoFChangeList< X > &dofs) override
 
template<typename... T, typename = std::enable_if_t<sizeof...(T) == N>>
void Interpolate (std::vector< X > &dofs, T... values) const
 Fill the degrees of freeedom by interpolating the specified constant function. More...
 
virtual std::size_t CodomainDimensionSize () const override
 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 override
 Populates an array with polynomial degree for each component of the result for the specified element. More...
 
void AddFunction (const std::shared_ptr< detail::DiscreteFunctionImpl< DIM, X, N >> &function)
 Adds the function to this function space. More...
 
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 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 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 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...
 

Static Protected Member Functions

static constexpr std::array< std::size_t, N > CreatePolyDegArray (std::size_t polynomialDegree)
 Create an array of size N filled with the specified polynomial degree. More...
 

Protected Attributes

PFEMesh< DIM > m_mesh
 Mesh the space is over. More...
 
FunctionList m_functions
 List of functions in the space. More...
 

Friends

template<std::size_t D, typename T , std::size_t M1, std::size_t... M>
PDiscreteFunctionSpace< D, T, detail::TemplateAdditionHack< M1, M... >::value > MakeDiscreteCartesianProductSpace (const PDiscreteFunctionSpace< D, T, M1 > &, const PDiscreteFunctionSpace< D, T, M > &...)
 

Related Functions

(Note that these are not member functions.)

template<std::size_t DIM, typename X , std::size_t M1, std::size_t... M>
PDiscreteFunctionSpace< DIM, X, detail::TemplateAdditionHack< M1, M... >::value > MakeDiscreteCartesianProductSpace (const PDiscreteFunctionSpace< DIM, X, M1 > &space, const PDiscreteFunctionSpace< DIM, X, M > &... spaces)
 Creates a Cartesian product space of the specified discrete function spaces. More...
 
template<std::size_t DIM, typename X , std::size_t N1, std::size_t N2>
PDiscreteFunctionSpace< DIM, X, N1+N2 > operator* (const PDiscreteFunctionSpace< DIM, X, N1 > &lhs, const PDiscreteFunctionSpace< DIM, X, N2 > &rhs)
 Creates a Cartesian product space of the specified discrete function spaces using the "multiplication" operator as the Cartesian product operator. More...
 

Detailed Description

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

Defines a space of discrete functions defined as the Cartesian product of discrete spaces.

Template Parameters
DIMThe dimension of the domain (and mesh)
XThe type of the result of the functions (double, std::complex, etc)
NThe size of the resulting vector space (N=1 for scalars)

Member Typedef Documentation

◆ CodomainType

template<std::size_t DIM, typename X , std::size_t N>
typedef X ptems::DiscreteFunctionSpace< DIM, X, N >::CodomainType
inherited

The type of the result of the functions.

◆ FunctionList

template<std::size_t DIM, typename X , std::size_t N>
typedef std::set<WeakFunction, std::owner_less<WeakFunction> > ptems::DiscreteFunctionSpace< DIM, X, N >::FunctionList
protectedinherited

List type (actually set) of weak pointers to functions constructed from this space.

◆ FunctionType

template<std::size_t DIM, typename X , std::size_t N>
typedef DiscreteFunction<DIM, X, N> ptems::DiscreteFunctionSpace< DIM, X, N >::FunctionType
inherited

The type of the functions in this space.

◆ WeakFunction

template<std::size_t DIM, typename X , std::size_t N>
typedef std::weak_ptr<detail::DiscreteFunctionImpl<DIM, X, N> > ptems::DiscreteFunctionSpace< DIM, X, N >::WeakFunction
protectedinherited

Weak pointer type to a function constructed from this space.

Member Function Documentation

◆ Adapt() [1/4]

template<std::size_t DIM, typename X , std::size_t N>
template<typename Container , typename = typename Container::value_type>
void ptems::DiscreteFunctionSpace< DIM, X, N >::Adapt ( const Container flags,
std::size_t  minimum = 1 
)
inlineinherited

Performs p-refinement on the space using a container of flags indicating whether to perform refinement (increment polynomial degree by one) or coarsen (decrement polynomial degree by one)

A container of flags can be one of the following:

  • A random access container of AdaptationType (i.e., supports operator[] and size()) containing an AdaptationType for each element
  • An associative container (map/unordered_map etc) of element number (std::size_t) to AdaptationType containing an AdaptationType for SOME elements, all other elements are left unrefined
  • A container of element numbers to perform refinement on, all other elements are left unrefined
Note
If the random access container size is less than the number of elements the remaining elements are left unrefined.
Attention
While the container of element numbers to refine can be ANY container type, it is strongly recommended to use a set-like container which can perform fast searching (e.g. std::set or std::unordered_set)
Template Parameters
ContainerContainer type containing AdaptationType or std::size_t, or map from std::size_t to AdaptationType
Parameters
flagsContainer of AdaptationType or std::size_t, or map from std::size_t to AdaptationType
minimumMinimum polynomial degree. Coarsen will not occur if the polynomial degree of ANY component will fall below this minimum

◆ Adapt() [2/4]

template<std::size_t DIM, typename X , std::size_t N>
template<typename Container >
void ptems::DiscreteFunctionSpace< DIM, X, N >::Adapt ( const Container refine,
const Container coarsen,
std::size_t  minimum = 1 
)
inlineinherited

Performs p-refinement on the space using a set of elements to refine (increment polynomial degree by one) and coarsen (decrement polynomial degree by one)

Attention
While the container of element numbers to refine/coarsen can be ANY container time, it is strongly recommended to use a set-like container which can perform fast searching (e.g. std::set or std::unordered_set)
Parameters
refineContainer of indices of elements to refine
coarsenContainer of indices of elements to coarsen
minimumMinimum polynomial degree. Coarsen will not occur if the polynomial degree of ANY component will fall below this minimum

◆ Adapt() [3/4]

template<std::size_t DIM, typename X , std::size_t N>
template<typename Func , typename = std::enable_if_t< std::is_same_v< decltype(std::declval<Func>()(std::declval<std::size_t>(), std::declval<std::array<std::size_t,N>&>())), bool>>
void ptems::DiscreteFunctionSpace< DIM, X, N >::Adapt ( Func  computePolydeg)
inlineinherited

Performs p-refinement on the space using a function which updates the polynomial degree of each element to the new polynomial degree.

Template Parameters
FuncA function type which takes a std::size_t (element number), and a REFERENCE to a std::array<std::size_t, N> containing the current polynomial degree of the function, and which should be updated with the new polynomial degree and returns a boolean indicating if the polynomial degree was changed; e.g., function of the form void operator()(std::size_t idx, std::array<std::size_t, N>& polydeg).
Parameters
computePolydegFunction which takes an element index in range [0, Mesh()->ElementCount()), and a C++ array of size N of polynominal degree of each element and updates the array with the new polynomial degree

◆ Adapt() [4/4]

template<std::size_t DIM, typename X , std::size_t N>
template<typename Func , typename = std::enable_if_t<std::is_same_v<decltype(std::declval<Func>()(std::declval<std::size_t>())), AdaptationType>>
void ptems::DiscreteFunctionSpace< DIM, X, N >::Adapt ( Func  flagFunction,
std::size_t  minimum = 1 
)
inlineinherited

Performs p-refinement on the space using a function for deciding which elements to refine.

Template Parameters
FuncA function type which takes a std::size_t, and returns an AdaptionType; e.g., of form AdaptionType operator()(std::size_t idx).
Parameters
flagFunctionFunction which takes an element index in range [0, Mesh()->ElementCount()), and returns an AdaptionType flag value indicating if the element polynomial degree should be refined (incremented by one), coarsened (decremented by one), or left 'as-is'
minimumMinimum polynomial degree. Coarsen will not occur if the polynomial degree of ANY component will fall below this minimum

◆ AddFunction()

template<std::size_t DIM, typename X , std::size_t N>
void ptems::DiscreteFunctionSpace< DIM, X, N >::AddFunction ( const std::shared_ptr< detail::DiscreteFunctionImpl< DIM, X, N >> &  function)
inlineprotectedinherited

Adds the function to this function space.

Parameters
functionFunction over this function space

◆ AttachMeshListener()

template<std::size_t DIM, typename X , std::size_t N>
void ptems::DiscreteFunctionSpace< DIM, X, N >::AttachMeshListener ( )
inlineinherited

Attach the listener for this space to the mesh.

Warning
This method MUST be called after construction. In general, functions spaces should not be constructed directly, but instead via MakeXXX functions which call this function automatically.

◆ BasisFunctions()

template<std::size_t DIM, typename X , std::size_t N>
DoFData<N, X> ptems::DiscreteFunctionSpace< DIM, X, N >::BasisFunctions ( std::size_t  elementIdx,
const Vector< DIM > &  pt 
)
inlineinherited

Gets the values of 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).

Note
For virtual element the returned values are the value projections of the basis functions
Parameters
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
The basis functions for each component of the space. The first index of the array is the component index and the second is the basis function index.

◆ CanIncrementPolynomialDegree() [1/2]

template<std::size_t DIM, typename X , std::size_t N>
bool ptems::DiscreteCartesianProductSpace< DIM, X, N >::CanIncrementPolynomialDegree ( ) const
inlineoverridevirtual

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 from ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ CanIncrementPolynomialDegree() [2/2]

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

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 , std::size_t N>
virtual std::size_t ptems::DiscreteFunctionSpace< DIM, X, N >::CodomainDimensionSize ( ) const
inlineoverrideprotectedvirtualinherited

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ ComputeBasisFunctions()

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

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ 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 
)
inlineprotectedvirtualinherited

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 , std::size_t N>
virtual std::vector<X>* ptems::DiscreteCartesianProductSpace< DIM, X, N >::ComputeDoFAtBasisProjection ( std::vector< X > *  dofs,
std::size_t  elementIdx 
)
inlineoverrideprotectedvirtual

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 from ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ ComputeDoFMapping()

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

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ ComputeGlobalDoFOffset()

template<std::size_t DIM, typename X , std::size_t N>
virtual void ptems::DiscreteCartesianProductSpace< DIM, X, N >::ComputeGlobalDoFOffset ( std::size_t dofOffset,
std::size_t  offset 
) const
inlineoverrideprotectedvirtual

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ ComputeGradBasisFunctions()

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

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ 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 
)
inlineprotectedvirtualinherited

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 
)
inlineprotectedvirtualinherited

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 
)
inlineprotectedvirtualinherited

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 , std::size_t N>
virtual void ptems::DiscreteCartesianProductSpace< DIM, X, N >::ComputeNumberLocalDoFs ( std::size_t  elementIdx,
std::size_t numberDofs 
) const
inlineoverrideprotectedvirtual

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ 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
inlineprotectedvirtualinherited

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

◆ CreatePolyDegArray()

template<std::size_t DIM, typename X , std::size_t N>
static constexpr std::array<std::size_t, N> ptems::DiscreteFunctionSpace< DIM, X, N >::CreatePolyDegArray ( std::size_t  polynomialDegree)
inlinestaticconstexprprotectedinherited

Create an array of size N filled with the specified polynomial degree.

Parameters
polynomialDegreeThe polynomial degree
Returns
Array filled with specified polynomial degree

◆ DoFMapping()

template<std::size_t DIM, typename X , std::size_t N>
DoFData<N, std::size_t> ptems::DiscreteFunctionSpace< DIM, X, N >::DoFMapping ( std::size_t  elementIdx) const
inlineinherited

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

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 local to global mapping of the degrees of freedom for
Returns
Mapping from local degree of freedom number to global degrees of freedom number for the specified element. The first index of the array is the component index and the second is the basis function index

◆ EstimateLocalAnalyticity()

template<std::size_t DIM, typename X , std::size_t N>
template<std::size_t M>
std::array<Vector<DIM>, N> ptems::DiscreteFunctionSpace< DIM, X, N >::EstimateLocalAnalyticity ( std::size_t  elementIdx,
const DiscreteFunction< DIM, X, M, N > &  func 
)
inlineinherited

Estimates the local analyticity of the element using Legendre coefficients.

Parameters
elementIdxElement to estimate analyticity
funcFunction to estimate analyticity of
Returns
Analyticity estimate, per variable and dimension, in range \((0,1]\), or negative if analyticity cannot be computed (not Legendre polynomials)

◆ EvaluateDoFAtBasis()

template<std::size_t DIM, typename X , std::size_t N>
DoFData<N, std::vector<X> > ptems::DiscreteFunctionSpace< DIM, X, N >::EvaluateDoFAtBasis ( std::size_t  elementIdx)
inlineinherited

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

Parameters
elementIdxThe element index to get the basis functions/dofs for
Returns
A vector for each degree of freedom containing the basis functions evaluated at that degree of freedom; e.g., for result result result(0,2)[1] returns the second basis of the first component evaluated at the third degree of freedom.

◆ EvaluateDoFAtBasisProjection()

template<std::size_t DIM, typename X , std::size_t N>
DoFData<N, std::vector<X> > ptems::DiscreteFunctionSpace< DIM, X, N >::EvaluateDoFAtBasisProjection ( std::size_t  elementIdx)
inlineinherited

Gets 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
elementIdxThe element index to get the basis functions/dofs for
Returns
A vector for each degree of freedom containing the value projection of the basis functions evaluated at that degree of freedom; e.g., for result result result(0,2)[1] returns the second basis of the first component evaluated at the third degree of freedom.

◆ FillPolynomialDegree()

template<std::size_t DIM, typename X , std::size_t N>
virtual std::size_t* ptems::DiscreteFunctionSpace< DIM, X, N >::FillPolynomialDegree ( std::size_t  element,
std::size_t nextPolydeg 
) const
inlineoverrideprotectedvirtualinherited

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()

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ Function()

template<std::size_t DIM, typename X , std::size_t N>
DiscreteFunction< DIM, X, N > ptems::DiscreteFunctionSpace< DIM, X, N >::Function
inherited

Creates a function in this function space.

Exceptions
std::logic_errorIf the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space.
Returns
Function over this function space

◆ GlobalDoFOffset()

template<std::size_t DIM, typename X , std::size_t N>
std::array<std::size_t,N+1> ptems::DiscreteFunctionSpace< DIM, X, N >::GlobalDoFOffset ( ) const
inlineinherited

Gets the offset into the global DoFs 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.
Returns
Array containing the offset into the global DoFs for the DoFs for each component (first N entries) and the total number of globals DoFs (last entry)

◆ GradBasisFunctions()

template<std::size_t DIM, typename X , std::size_t N>
DoFData<N, FuncAndGradData<DIM,X> > ptems::DiscreteFunctionSpace< DIM, X, N >::GradBasisFunctions ( std::size_t  elementIdx,
const Vector< DIM > &  pt 
)
inlineinherited

Gets the derivative of the basis functions in each coordinate direction for each component at the specified point in 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).

Note
For virtual element the returned values are the value and gradient projections of the basis functions
Parameters
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
The derivatives of the basis functions for each component of the space. The first index of the array is the component index, the second is the basis function index, and the third is the index of the derivative direction.

◆ GradValueBasisFunctions()

template<std::size_t DIM, typename X , std::size_t N>
DoFData<N, FuncAndGradData<DIM,X> > ptems::DiscreteFunctionSpace< DIM, X, N >::GradValueBasisFunctions ( std::size_t  elementIdx,
const Vector< DIM > &  pt 
)
inlineinherited

Gets the derivative of the value of basis functions in each coordinate direction for each component at the specified point in 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).

Note
For virtual element the returned values are the value projection and gradient of the value projection. This is different to GradBasisFunctions
Parameters
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
The derivatives of the basis functions for each component of the space. The first index of the array is the component index, the second is the basis function index, and the third is the index of the derivative direction.

◆ HessianBasisFunctions()

template<std::size_t DIM, typename X , std::size_t N>
DoFData<N, FuncGradAndHessianData<DIM,X> > ptems::DiscreteFunctionSpace< DIM, X, N >::HessianBasisFunctions ( std::size_t  elementIdx,
const Vector< DIM > &  pt 
)
inlineinherited

Gets the second order derivatives of the basis functions for each component at the specified point in 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).

Exceptions
std::logic_errorIf second order derivatives are not supported by the space
Parameters
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
The second order derivatives of the basis functions for each component of the space. The first index of the array is the component index, the second is the basis function index, and the third/fourth are the indices of the derivatives

◆ Interpolant()

template<std::size_t DIM, typename X , std::size_t N>
template<typename... T, typename >
DiscreteFunction< DIM, X, N > ptems::DiscreteFunctionSpace< DIM, X, N >::Interpolant ( T...  values)
inherited

Creates a function in this function space, with initial constant function value.

Template Parameters
TThe type of the constant values, must be convertible to X. Must have 1 or N arguments
Parameters
valuesConstant value for all components of the function (if one argument specified), or the individual values for each component of the function (if N arguments specified)

◆ Interpolate() [1/2]

template<std::size_t DIM, typename X , std::size_t N>
template<typename... T, typename = std::enable_if_t<sizeof...(T) == N>>
void ptems::DiscreteFunctionSpace< DIM, X, N >::Interpolate ( std::vector< X > &  dofs,
T...  values 
) const
inlineprotectedinherited

Fill the degrees of freeedom by interpolating the specified constant function.

Template Parameters
TThe type of the constants
Parameters
dofsThe degrees of freedom to fill
valuesThe constants (one constant per component of the space)

◆ Interpolate() [2/2]

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

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().

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ IsContinuous()

template<std::size_t DIM, typename X , std::size_t N>
bool ptems::DiscreteCartesianProductSpace< DIM, X, N >::IsContinuous ( ) const
inlineoverridevirtual

Gets if the functions space is continuous over element boundaries.

Returns
true If continuous, false otherwise

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ Mesh()

template<std::size_t DIM, typename X , std::size_t N>
const PFEMesh<DIM>& ptems::DiscreteFunctionSpace< DIM, X, N >::Mesh ( ) const
inlineinherited

Gets the mesh the function space is over.

Returns
The mesh

◆ MeshChanged() [1/2]

template<std::size_t DIM, typename X , std::size_t N>
void ptems::DiscreteFunctionSpace< DIM, X, N >::MeshChanged ( const PFEMesh< DIM > &  previousMesh,
const typename FEMesh< DIM >::Modifications &  changes 
)
inlineoverrideinherited

Event notification from mesh that the mesh has been changed.

Exceptions
std::logic_errorIf the new mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space.
Attention
This method is used as a callback from FEMesh; hence, it should generally not be called in client code.
Parameters
previousMeshThe mesh before refinement
changesA FEMesh<DIM>::Modifications structure defining the change to the mesh

◆ MeshChanged() [2/2]

template<std::size_t DIM>
virtual void ptems::FEMesh< DIM >::MeshChangeListener::MeshChanged ( const std::shared_ptr< FEMesh< DIM >> &  previousMesh,
const typename FEMesh< DIM >::Modifications changes 
)
pure virtualinherited

Event notification from mesh that the mesh has been changed.

Exceptions
std::logic_errorIf the new mesh is invalid
Parameters
previousMeshThe mesh before refinement
changesA FEMesh<DIM>::Modifications structure defining the change to the mesh

◆ NumberGlobalDoFs()

template<std::size_t DIM, typename X , std::size_t N>
virtual std::size_t ptems::DiscreteCartesianProductSpace< DIM, X, N >::NumberGlobalDoFs ( ) const
inlineoverridevirtual

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ NumberLocalDoFs()

template<std::size_t DIM, typename X , std::size_t N>
std::array<std::size_t,N> ptems::DiscreteFunctionSpace< DIM, X, N >::NumberLocalDoFs ( std::size_t  elementIdx) const
inlineinherited

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
Returns
Local number of degrees of freedom for the specified element per component

◆ OnMeshChanged()

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

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

Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ PolynomialDegree() [1/4]

template<std::size_t DIM, typename X , std::size_t N>
const std::array<std::size_t, N>& ptems::DiscreteCartesianProductSpace< DIM, X, N >::PolynomialDegree ( std::size_t  element) const
inlineoverridevirtual

Gets the polynomial degree of the function space for the specified element.

Parameters
elementThe element to get the polynomial space from
Returns
The polynomial degree of the function space in each component

Implements ptems::DiscreteFunctionSpace< DIM, X, N >.

◆ PolynomialDegree() [2/4]

template<std::size_t DIM, typename X , std::size_t N>
virtual void ptems::DiscreteFunctionSpace< DIM, X, N >::PolynomialDegree ( std::size_t  element,
const std::array< std::size_t, N > &  polydeg 
)
inlinevirtualinherited

Sets the polynomial degree of the function space for the specified element.

May be ignored if SupportsVariablePolynomialDegree() is false

Parameters
elementThe element to get the polynomial space from
polydegThe polynomial degree of the function space in each component

◆ PolynomialDegree() [3/4]

template<std::size_t DIM, typename X , std::size_t N>
void ptems::DiscreteFunctionSpace< DIM, X, N >::PolynomialDegree ( std::size_t  element,
std::size_t  polydeg 
)
inlineinherited

Sets the polynomial degree of the function space for the specified element.

May be ignored if SupportsvariablePolynomialDegree() is false

Parameters
elementThe element to get the polynomial space from
polydegThe polynomial degree of the function space for all components

◆ PolynomialDegree() [4/4]

template<std::size_t DIM, typename X , std::size_t N>
template<typename... T, typename = std::enable_if_t<sizeof...(T) == N && (sizeof...(T) > 1)>>
void ptems::DiscreteFunctionSpace< DIM, X, N >::PolynomialDegree ( std::size_t  element,
T...  polydeg 
)
inlineinherited

Sets the polynomial degree of the function space for the specified element.

May be ignored if SupportsvariablePolynomialDegree() is false

Parameters
elementThe element to get the polynomial space from
polydegThe polynomial degree of the function space in each component

◆ PolynomialDegreeFunction()

template<std::size_t DIM, typename X , std::size_t N>
SpacePolynomialDegree ptems::DiscreteFunctionSpace< DIM, X, N >::PolynomialDegreeFunction ( std::size_t  component = 0) const
inlineinherited

Gets the polynomial degree of the function space for the specified component.

Parameters
componentThe component to get the polynomial space from
Returns
The polynomial degree of the function space in each component

◆ 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 
)
inlineprotectedvirtualinherited

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 
)
inlineprotectedvirtualinherited

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

◆ StronglyImposedDoFs()

template<std::size_t DIM, typename X , std::size_t N>
std::unordered_map<std::size_t, X> ptems::DiscreteFunctionSpace< DIM, X, N >::StronglyImposedDoFs ( ) const
inlineinherited

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

Returns
dofs Set of global dof number and value for strongly imposed Dirichlet BC

◆ SupportsVariablePolynomialDegree() [1/2]

template<std::size_t DIM, typename X , std::size_t N>
bool ptems::DiscreteCartesianProductSpace< DIM, X, N >::SupportsVariablePolynomialDegree ( ) const
inlineoverridevirtual

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 from ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ SupportsVariablePolynomialDegree() [2/2]

template<std::size_t DIM, typename X , std::size_t N>
bool ptems::DiscreteCartesianProductSpace< DIM, X, N >::SupportsVariablePolynomialDegree ( [[maybe_unused] ] std::size_t  component) const
inlineoverridevirtual

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 from ptems::DiscreteFunctionSpaceInterface< DIM, X >.

◆ UniformPolynomialDegree() [1/3]

template<std::size_t DIM, typename X , std::size_t N>
virtual void ptems::DiscreteFunctionSpace< DIM, X, N >::UniformPolynomialDegree ( const std::array< std::size_t, N > &  polydeg)
inlinevirtualinherited

Sets the polynomial degree of the function space for all elements.

May be ignored if SupportsvariablePolynomialDegree() is false

Parameters
polydegThe polynomial degree of the function space in each component

◆ UniformPolynomialDegree() [2/3]

template<std::size_t DIM, typename X , std::size_t N>
void ptems::DiscreteFunctionSpace< DIM, X, N >::UniformPolynomialDegree ( std::size_t  polydeg)
inlineinherited

Sets the polynomial degree of the function space for all elements.

May be ignored if SupportsvariablePolynomialDegree() is false

Parameters
polydegThe polynomial degree of the function space for all components

◆ UniformPolynomialDegree() [3/3]

template<std::size_t DIM, typename X , std::size_t N>
template<typename... T, typename = std::enable_if_t<sizeof...(T) == N && (sizeof...(T) > 1)>>
void ptems::DiscreteFunctionSpace< DIM, X, N >::UniformPolynomialDegree ( T...  polydeg)
inlineinherited

Sets the polynomial degree of the function space for all elements.

May be ignored if SupportsvariablePolynomialDegree() is false

Parameters
polydegThe polynomial degree of the function space in each component

Friends And Related Function Documentation

◆ MakeDiscreteCartesianProductSpace()

template<std::size_t DIM, typename X , std::size_t M1, std::size_t... M>
PDiscreteFunctionSpace< DIM, X, detail::TemplateAdditionHack< M1, M... >::value > MakeDiscreteCartesianProductSpace ( const PDiscreteFunctionSpace< DIM, X, M1 > &  space,
const PDiscreteFunctionSpace< DIM, X, M > &...  spaces 
)
related

Creates a Cartesian product space of the specified discrete function spaces.

Exceptions
std::invalid_argumentIf any of the spaces do not use the same mesh
Template Parameters
M1The size of the resulting vector space of the first space to Cartesian product
MThe sizes of the resulting vector space of the second (and subsequent) spaces to Cartesian product
Parameters
spaceThe first space to Cartesian product
spacesThe second (and subsequent) spaces to Cartesian product

◆ operator*()

template<std::size_t DIM, typename X , std::size_t N1, std::size_t N2>
PDiscreteFunctionSpace< DIM, X, N1+N2 > operator* ( const PDiscreteFunctionSpace< DIM, X, N1 > &  lhs,
const PDiscreteFunctionSpace< DIM, X, N2 > &  rhs 
)
related

Creates a Cartesian product space of the specified discrete function spaces using the "multiplication" operator as the Cartesian product operator.

Exceptions
std::invalid_argumentIf spaces do not use the same mesh
Template Parameters
N1The size of the resulting vector space of the first space to Cartesian product
N2The sizes of the resulting vector space of the second spaces to Cartesian product
Parameters
lhsThe first space to Cartesian product
rhsThe second space to Cartesian product

Member Data Documentation

◆ CodomainDimension

template<std::size_t DIM, typename X , std::size_t N>
constexpr std::size_t ptems::DiscreteFunctionSpace< DIM, X, N >::CodomainDimension = N
staticconstexprinherited

The size of the resulting vector space for the functions.

◆ DomainDimension

template<std::size_t DIM, typename X , std::size_t N>
constexpr std::size_t ptems::DiscreteFunctionSpace< DIM, X, N >::DomainDimension = DIM
staticconstexprinherited

The dimension of the spacial domain the functions are over.

◆ m_functions

template<std::size_t DIM, typename X , std::size_t N>
FunctionList ptems::DiscreteFunctionSpace< DIM, X, N >::m_functions
protectedinherited

List of functions in the space.

◆ m_mesh

template<std::size_t DIM, typename X , std::size_t N>
PFEMesh<DIM> ptems::DiscreteFunctionSpace< DIM, X, N >::m_mesh
protectedinherited

Mesh the space is over.


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