Defines a space of discrete functions defined as the Cartesian product of discrete spaces. More...
#include <ptems/space.hpp>
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_t > | DoFMapping (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... | |
T | shared_from_this (T... args) |
T | 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_t * | ComputeDoFMapping (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_t * | FillPolynomialDegree (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... | |
Defines a space of discrete functions defined as the Cartesian product of discrete spaces.
DIM | The dimension of the domain (and mesh) |
X | The type of the result of the functions (double, std::complex, etc) |
N | The size of the resulting vector space (N=1 for scalars) |
|
inherited |
The type of the result of the functions.
|
protectedinherited |
List type (actually set) of weak pointers to functions constructed from this space.
|
inherited |
The type of the functions in this space.
|
protectedinherited |
Weak pointer type to a function constructed from this space.
|
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:
Container | Container type containing AdaptationType or std::size_t, or map from std::size_t to AdaptationType |
flags | Container of AdaptationType or std::size_t, or map from std::size_t to AdaptationType |
minimum | Minimum polynomial degree. Coarsen will not occur if the polynomial degree of ANY component will fall below this minimum |
|
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)
|
inlineinherited |
Performs p-refinement on the space using a function which updates the polynomial degree of each element to the new polynomial degree.
Func | A 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) . |
computePolydeg | Function 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 |
|
inlineinherited |
Performs p-refinement on the space using a function for deciding which elements to refine.
Func | A function type which takes a std::size_t, and returns an AdaptionType; e.g., of form AdaptionType operator()(std::size_t idx) . |
flagFunction | Function 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' |
minimum | Minimum polynomial degree. Coarsen will not occur if the polynomial degree of ANY component will fall below this minimum |
|
inlineprotectedinherited |
Adds the function to this function space.
function | Function over this function space |
|
inlineinherited |
Attach the listener for this space to the mesh.
|
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).
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
|
inlineoverridevirtual |
Gets whether the space can be modified by incrementing/decrementing the polynomial degree uniformly for any component.
Reimplemented from ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlinevirtualinherited |
Gets whether the space can be modified by incrementing/decrementing the polynomial degree uniformly for the specified codomain component.
component | Index of component to check |
|
inlineoverrideprotectedvirtualinherited |
Gets ths size of the codomain of the fuctions in the space.
Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
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
basis | Pointer 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. |
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
basis+x
, where x
is the sum of the array entries computed by ComputeNumberLocalDoFs Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlineprotectedvirtualinherited |
Computes the basis functions of the space evaluated by the functional for each degree of freedom.
dofs | Pointer 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 |
elementIdx | The element index to get the basis functions/dofs for |
dofs+x
, where x
is the sum of the array entries computed by ComputeNumberLocalDoFs
|
inlineoverrideprotectedvirtual |
Computes the value projection of the basis functions of the space evaluated by the functional for each degree of freedom.
dofs | Pointer 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 |
elementIdx | The element index to get the basis functions/dofs for |
dofs+x
, where x
is the sum of the array entries computed by ComputeNumberLocalDoFs Reimplemented from ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
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
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
std::out_of_range | if elementIdx is not less than the number of elements in the underlying mesh |
dofMapping | Pointer 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. |
elementIdx | Index of the element to get the local to global mapping of the degrees of freedom for |
offset | Offset to add to the global DoF numbers |
dofMapping+x
, where x
is the sum of the array entries computed by ComputeNumberLocalDoFs Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlineoverrideprotectedvirtual |
Computes the offset onto the global degrees of freedom for each component.
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
dofOffset | Pointer 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) |
offset | Offset to add to the global DoF numbers |
Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
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
basis | Pointer 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. |
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
basis+x
, where x
is the sum of the array entries computed by ComputeNumberLocalDoFs Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
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
basis | Pointer 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. |
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
basis+x
, where x
is the sum of the array entries computed by ComputeNumberLocalDoFs
|
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
std::logic_error | If second order derivatives are not supported by the space |
basis | Pointer 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. |
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
basis+x
, where x
is the sum of the array entries computed by ComputeNumberLocalDoFs
|
inlineprotectedvirtualinherited |
Computes the local analyticity estimate of the element using Legendre coefficients.
elementIdx | Element to estimate analyticity |
analyticity | Pointer to the array to store analyticity estimate, per variable, in range \((0,1]\) (or negative if analyticity cannot be computed) |
dofs | Pointer to array containing degrees of freedom |
|
inlineoverrideprotectedvirtual |
Gets the number of degrees of freedom for an element in the space per component.
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
std::out_of_range | if elementIdx is not less than the number of elements in the underlying mesh |
elementIdx | Index of the element to get the number of local degrees of freedom for |
numberDofs | Pointer 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 >.
|
inlineprotectedvirtualinherited |
Computes the set of global DoFs with strongly imposed Dirichlet boundary conditions.
dofs | Set to update with global dof number and value for strongly imposed Dirichlet BC |
offset | Offset to add to the global DoF numbers |
|
inlinestaticconstexprprotectedinherited |
Create an array of size N filled with the specified polynomial degree.
polynomialDegree | The polynomial degree |
|
inlineinherited |
Gets the mapping between local dof numbers for an element and the global dof numbers.
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
std::out_of_range | if elementIdx is not less than the number of elements in the underlying mesh |
elementIdx | Index of the element to get the local to global mapping of the degrees of freedom for |
|
inlineinherited |
Estimates the local analyticity of the element using Legendre coefficients.
elementIdx | Element to estimate analyticity |
func | Function to estimate analyticity of |
|
inlineinherited |
Gets the basis functions of the space evaluated by the functional for each degree of freedom.
elementIdx | The element index to get the basis functions/dofs for |
result
result(0,2)[1]
returns the second basis of the first component evaluated at the third degree of freedom.
|
inlineinherited |
Gets the value projection of the basis functions of the space evaluated by the functional for each degree of freedom.
elementIdx | The element index to get the basis functions/dofs for |
result
result(0,2)[1]
returns the second basis of the first component evaluated at the third degree of freedom.
|
inlineoverrideprotectedvirtualinherited |
Populates an array with polynomial degree for each component of the result for the specified element.
element | The index of the element to get the polynomial degree for |
nextPolydeg | Point into an array to store the polynomial degree. The array will have enough space to store CodomainDimensionSize() values |
Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inherited |
Creates a function in this function space.
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
|
inlineinherited |
Gets the offset into the global DoFs for each component.
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
|
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).
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
|
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).
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
|
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).
std::logic_error | If second order derivatives are not supported by the space |
elementIdx | The element index to get the basis functions for |
pt | The local point in the reference element or bounding box to get the basis functions for |
|
inherited |
Creates a function in this function space, with initial constant function value.
T | The type of the constant values, must be convertible to X. Must have 1 or N arguments |
values | Constant 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) |
|
inlineprotectedinherited |
Fill the degrees of freeedom by interpolating the specified constant function.
T | The type of the constants |
dofs | The degrees of freedom to fill |
values | The constants (one constant per component of the space) |
|
inlineoverrideprotectedvirtual |
Fill the degrees of freedom for this space by interpolating the specified constant function.
std::logic_error | If dofsEnd-dofsBegin is less than the number of degrees of freedom to be populated |
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()
dofsBegin | Iterator pointing to the start of the degrees of freedom to fill |
dofsEnd | Iterator pointing to the end of the degrees of freedom to fill |
value | Pointer to array containing the constants for each component |
dofsBegin+NumberGlobalDoFs()
. Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlineoverridevirtual |
Gets if the functions space is continuous over element boundaries.
Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlineinherited |
Gets the mesh the function space is over.
|
inlineoverrideinherited |
Event notification from mesh that the mesh has been changed.
std::logic_error | If the new mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
previousMesh | The mesh before refinement |
changes | A FEMesh<DIM>::Modifications structure defining the change to the mesh |
|
pure virtualinherited |
Event notification from mesh that the mesh has been changed.
std::logic_error | If the new mesh is invalid |
previousMesh | The mesh before refinement |
changes | A FEMesh<DIM>::Modifications structure defining the change to the mesh |
|
inlineoverridevirtual |
Gets the global number of degrees of freedom for the space.
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
Implements ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlineinherited |
Gets the number of degrees of freedom for an element in the space per component.
std::logic_error | If the underlying mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
std::out_of_range | if elementIdx is not less than the number of elements in the underlying mesh |
elementIdx | Index of the element to get the number of local degrees of freedom for |
|
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.
std::logic_error | If the new mesh is invalid for the type of space; e.g., nonconforming mesh for a conforming space. |
previousMesh | The mesh before refinement |
changes | A FEMesh<DIM>::Modifications structure defining the change to the mesh |
dofs | Map (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 >.
|
inlineoverridevirtual |
Gets the polynomial degree of the function space for the specified element.
element | The element to get the polynomial space from |
Implements ptems::DiscreteFunctionSpace< DIM, X, N >.
|
inlinevirtualinherited |
Sets the polynomial degree of the function space for the specified element.
May be ignored if SupportsVariablePolynomialDegree() is false
element | The element to get the polynomial space from |
polydeg | The polynomial degree of the function space in each component |
|
inlineinherited |
Sets the polynomial degree of the function space for the specified element.
May be ignored if SupportsvariablePolynomialDegree() is false
element | The element to get the polynomial space from |
polydeg | The polynomial degree of the function space for all components |
|
inlineinherited |
Sets the polynomial degree of the function space for the specified element.
May be ignored if SupportsvariablePolynomialDegree() is false
element | The element to get the polynomial space from |
polydeg | The polynomial degree of the function space in each component |
|
inlineinherited |
Gets the polynomial degree of the function space for the specified component.
component | The component to get the polynomial space from |
|
inlineprotectedvirtualinherited |
Set the polynomial degree of the specified elements and update the DoFs of the specified functions.
elementPolydegs | (Ordered) map of element number to array (of length N) of new polynomial degrees for each component |
dofs | Map (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 |
|
inlineprotectedvirtualinherited |
Set the polynomial degree of the all elements and update the DoFs of the specified functions.
polydeg | Array (of length N) of new polynomial degrees for each component |
dofs | Map (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 |
|
inlineinherited |
Computes the set of global DoFs with strongly imposed Dirichlet boundary conditions.
|
inlineoverridevirtual |
Gets whether the space can be modified with varying polynomial degree for each element for any component.
Reimplemented from ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlineoverridevirtual |
Gets whether the space can be modified with varying polynomial degree for each element for the specified codomain component.
component | Index of component to check |
Reimplemented from ptems::DiscreteFunctionSpaceInterface< DIM, X >.
|
inlinevirtualinherited |
Sets the polynomial degree of the function space for all elements.
May be ignored if SupportsvariablePolynomialDegree() is false
polydeg | The polynomial degree of the function space in each component |
|
inlineinherited |
Sets the polynomial degree of the function space for all elements.
May be ignored if SupportsvariablePolynomialDegree() is false
polydeg | The polynomial degree of the function space for all components |
|
inlineinherited |
Sets the polynomial degree of the function space for all elements.
May be ignored if SupportsvariablePolynomialDegree() is false
polydeg | The polynomial degree of the function space in each component |
|
related |
Creates a Cartesian product space of the specified discrete function spaces.
std::invalid_argument | If any of the spaces do not use the same mesh |
M1 | The size of the resulting vector space of the first space to Cartesian product |
M | The sizes of the resulting vector space of the second (and subsequent) spaces to Cartesian product |
space | The first space to Cartesian product |
spaces | The second (and subsequent) spaces to Cartesian product |
|
related |
Creates a Cartesian product space of the specified discrete function spaces using the "multiplication" operator as the Cartesian product operator.
std::invalid_argument | If spaces do not use the same mesh |
N1 | The size of the resulting vector space of the first space to Cartesian product |
N2 | The sizes of the resulting vector space of the second spaces to Cartesian product |
lhs | The first space to Cartesian product |
rhs | The second space to Cartesian product |
|
staticconstexprinherited |
The size of the resulting vector space for the functions.
|
staticconstexprinherited |
The dimension of the spacial domain the functions are over.
|
protectedinherited |
List of functions in the space.
|
protectedinherited |
Mesh the space is over.