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

Class for performing \(hp\)-adaptive refinement, where the selection on whether to perform \(h\)- or \(p\)-refinement is based on the decay of the Legendre coefficients of the numerical solution; cf. More...

#include <ptems/hpstrategies.hpp>

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

Public Member Functions

 HPLegendreDecayStrategy (const PDiscreteFunctionSpace< DIM, X, N > &space, double threshold=0.25, std::size_t minP=1)
 Creates an object for performing \(hp\)-adaptive refinement based on the decay of Legendre coefficients. More...
 
const PDiscreteFunctionSpace< DIM, double, N > & Space ()
 Gets the finite element space the refinement strategy is on. More...
 
FEMesh< DIM >::Modifications AdaptWithStrategyAndDefaultSmoothing (const Estimators &estimators, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor, and with default mesh smoothing options (default constructed MeshSmoother). More...
 
FEMesh< DIM >::Modifications AdaptWithStrategyAndDefaultSmoothing (const Estimators &estimators, const DiscreteFunction< DIM, double, N, N > &func, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified function defined on the space, and with default mesh smoothing options (default constructed MeshSmoother). More...
 
FEMesh< DIM >::Modifications AdaptWithStrategyAndDefaultSmoothing (const Estimators &estimators, const DiscreteFuncs &funcs, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified functions defined on the space, and with default mesh smoothing options (default constructed MeshSmoother). More...
 
FEMesh< DIM >::Modifications AdaptConformingWithStrategy (const Estimators &estimators, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor, with mesh smoothing options for conforming meshes (i.e., MeshSmoother::Conforming()). More...
 
FEMesh< DIM >::Modifications AdaptConformingWithStrategy (const Estimators &estimators, const DiscreteFunction< DIM, double, N, N > &func, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified function defined on the space with mesh smoothing options for conforming meshes (i.e., MeshSmoother::Conforming()). More...
 
FEMesh< DIM >::Modifications AdaptConformingWithStrategy (const Estimators &estimators, const DiscreteFuncs &funcs, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified functions defined on the space with mesh smoothing options for conforming meshes (i.e., MeshSmoother::Conforming()). More...
 
FEMesh< DIM >::Modifications AdaptWithStrategy (const Estimators &estimators, MeshSmoother smoother, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor and with the specified mesh smoothing options. More...
 
FEMesh< DIM >::Modifications AdaptWithStrategy (const Estimators &estimators, const DiscreteFunction< DIM, double, N, N > &func, MeshSmoother smoother, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified function defined on the space and with the specified mesh smoothing options. More...
 
FEMesh< DIM >::Modifications AdaptWithStrategy (const Estimators &estimators, const DiscreteFuncs &funcs, MeshSmoother smoother, Args... args)
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified functions defined on the space and with the specified mesh smoothing options. More...
 
FEMesh< DIM >::Modifications Adapt (const Container &flags, const Estimators &estimators, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and using a container of flags indicating the refinement to perform. More...
 
FEMesh< DIM >::Modifications Adapt (const Container &flags, const Estimators &estimators, const DiscreteFunction< DIM, double, N, N > &func, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and using a container of flags indicating the refinement to perform. More...
 
FEMesh< DIM >::Modifications Adapt (const Container &flags, const Estimators &estimators, const DiscreteFuncs &funcs, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and using a container of flags indicating the refinement to perform. More...
 
FEMesh< DIM >::Modifications Adapt (const Container &refine, const Container &coarsen, const Estimators &estimators, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a set of elements to refine and coarsen. More...
 
FEMesh< DIM >::Modifications Adapt (const Container &refine, const Container &coarsen, const Estimators &estimators, const DiscreteFunction< DIM, double, N, N > &func, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a set of elements to refine and coarsen. More...
 
FEMesh< DIM >::Modifications Adapt (const Container &refine, const Container &coarsen, const Estimators &estimators, const DiscreteFuncs &funcs, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a set of elements to refine and coarsen. More...
 
FEMesh< DIM >::Modifications Adapt (Func flagFunction, const Estimators &estimators, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a function for deciding which elements to refine. More...
 
FEMesh< DIM >::Modifications Adapt (Func flagFunction, const Estimators &estimators, const DiscreteFunction< DIM, double, N, N > &func, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a function for deciding which elements to refine. More...
 
FEMesh< DIM >::Modifications Adapt (Func flagFunction, const Estimators &estimators, const DiscreteFuncs &funcs, MeshSmoother smoother=MeshSmoother{})
 Perform \(hp\)-adaptive refinement on the contained space and mesh using a function for deciding which elements to refine. More...
 

Protected Member Functions

virtual bool ShouldPolydegAdapt (const HPRefinementData< DIM, X, N > &data, std::size_t element, bool refine) override
 Method called by the adaptation routines to decide if \(p\)-refinement should be performed for an element. More...
 
virtual bool ShouldPolydegAdapt (const HPRefinementData< DIM, double, N > &data, std::size_t element, bool refine)=0
 Method called by the adaptation routines to decide if \(p\)-refinement should be performed for an element. More...
 
virtual void PostProcessSpaceRefine ([[maybe_unused]] const HPRefinementData< DIM, double, N > &data,[[maybe_unused]] const std::vector< AdaptationType > &pRefinements)
 Method called by the adaptation routines after \(p\)-refinement, but before \(h\)-refinement. More...
 
virtual void PostProcessMeshRefine ([[maybe_unused]] const HPRefinementData< DIM, double, N > &data,[[maybe_unused]] const typename FEMesh< DIM >::Modifications &changes,[[maybe_unused]] const std::vector< AdaptationType > &pRefinements)
 Method called by the adaptation routines after both \(h\)-refinement and \(p\)-refinement. More...
 

Detailed Description

template<std::size_t DIM, typename X = double, std::size_t N = 1>
class ptems::HPLegendreDecayStrategy< DIM, X, N >

Class for performing \(hp\)-adaptive refinement, where the selection on whether to perform \(h\)- or \(p\)-refinement is based on the decay of the Legendre coefficients of the numerical solution; cf.

[4].

For each element \(K\in\mathcal{T}_h\) marked for refinement the algorithm performs \(p\)-refinement if \(\theta < \gamma\), where \(\theta>0\) is an estimate of the smoothness, with \(\theta\) tending to zero for smooth functions, and \(\gamma\) is a user specified threshold.

Warning
If the polynomials used by the space are not Legendre polynomials this method will NEVER perform \(p\)-refinement.

Usage: An instance of HPLegendreDecayStrategy should be constructed when the initial mesh and space are constructed, and then the various Adapt* functions should be called to perform adaptation on BOTH the mesh and space:

std::size_t polydeg = 2;
auto Vh = ptems::MakeDGSpace(mesh, polydeg);
auto hpRefinement = ptems::CreateHPLegendreDecayStrategy(Vh, 0.25, polydeg);
for(std::size_t i = 0; i != 20; ++i)
{
// ... Perform solve
// ... Compute error indicators:
std::vector<double> indicators(mesh->ElementCount());
// Refine using a marking strategy functor:
hpRefinement->AdaptWithStrategy<ptems::FixedFractionMarking>(indicators, uh, 0.25);
// Refine using custom functor (or similar for container of flags)
hpRefinement->Adapt([](std::size_t element)
{
// return value from AdaptationType enum
}, indicators, uh);
}
static DomainDefinition2D Square()
Create a unit square computational domain.
Definition: domaindefinition.hpp:1692
PFEMesh< 1 > CreateUniformMesh(const ParallelComm &mpi, double min, double max, std::size_t numElements)
Create a 1 dimensional finite element mesh by uniform subdivision of the specified interval.
PHPRefinementStrategy< DIM, X, N > CreateHPLegendreDecayStrategy(const PDiscreteFunctionSpace< DIM, X, N > &space, double threshold=0.25, std::size_t minP=1)
Creates an object for performing -adaptive refinement based on the estimated smoothness of the numeri...
Definition: hpstrategies.hpp:1174
PDiscreteFunctionSpace< DIM, X, N > MakeDGSpace(const PFEMesh< DIM > &mesh, std::size_t polynomialDegree, TensorProductPolynomials tensorProductPolynomials=TensorProductPolynomials::HyperplaneExtrudedElementsOnly)
Create a discontinuous Galerkin finite element space over the specified mesh with specified polynomia...
Definition: discontinuousfem.hpp:596
Functor object that can be passed to FEMesh::Adapt(Func, bool) to perform Fixed Fraction marking for ...
Definition: refinements.hpp:448
Note
It is simpler to call CreateHPLegendreDecayStrategy as this automatically deduces template parameters.
Template Parameters
DIMThe dimension of the mesh and spaces to refine
XThe type of the result of the functions (double, std::complex, etc) for each space
NThe size of the resulting vector space (N=1 for scalars) for each space

Constructor & Destructor Documentation

◆ HPLegendreDecayStrategy()

template<std::size_t DIM, typename X = double, std::size_t N = 1>
ptems::HPLegendreDecayStrategy< DIM, X, N >::HPLegendreDecayStrategy ( const PDiscreteFunctionSpace< DIM, X, N > &  space,
double  threshold = 0.25,
std::size_t  minP = 1 
)
inline

Creates an object for performing \(hp\)-adaptive refinement based on the decay of Legendre coefficients.

Parameters
spaceThe function space to refine (containing the mesh to refine)
thresholdThreshold \(\gamma\) specifying the value the decay should be below to perform \(p\)-refinement (default=0.25)
minPMinimum polynomial degree. Coarsen will not occur if the polynomial degree of ANY component will fall below this minimum (default = 1)

Member Function Documentation

◆ Adapt() [1/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( const Container flags,
const Estimators &  estimators,
const DiscreteFuncs &  funcs,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and using a container of flags indicating the refinement to perform.

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
Template Parameters
ContainerContainer type containing AdaptationType or std::size_t, or map from std::size_t to AdaptationType
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
DiscreteFuncsType of the list of functions, must be random access container.
Parameters
flagsContainer of AdaptationType or std::size_t, or map from std::size_t to AdaptationType
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
funcsList of functions defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [2/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( const Container flags,
const Estimators &  estimators,
const DiscreteFunction< DIM, X, N, N > &  func,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and using a container of flags indicating the refinement to perform.

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
Template Parameters
ContainerContainer type containing AdaptationType or std::size_t, or map from std::size_t to AdaptationType
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
Parameters
flagsContainer of AdaptationType or std::size_t, or map from std::size_t to AdaptationType
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
funcFunction defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [3/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( const Container flags,
const Estimators &  estimators,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and using a container of flags indicating the refinement to perform.

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
Template Parameters
ContainerContainer type containing AdaptationType or std::size_t, or map from std::size_t to AdaptationType
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
Parameters
flagsContainer of AdaptationType or std::size_t, or map from std::size_t to AdaptationType
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [4/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( const Container refine,
const Container coarsen,
const Estimators &  estimators,
const DiscreteFuncs &  funcs,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a set of elements to refine and coarsen.

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)
Template Parameters
ContainerContainer type for refinement and coarsen flags
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
DiscreteFuncsType of the list of functions, must be random access container.
Parameters
refineContainer of indices of elements to refine
coarsenContainer of indices of elements to coarsen
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
funcsList of functions defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [5/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( const Container refine,
const Container coarsen,
const Estimators &  estimators,
const DiscreteFunction< DIM, X, N, N > &  func,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a set of elements to refine and coarsen.

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)
Template Parameters
ContainerContainer type for refinement and coarsen flags
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
Parameters
refineContainer of indices of elements to refine
coarsenContainer of indices of elements to coarsen
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
funcFunction defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [6/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( const Container refine,
const Container coarsen,
const Estimators &  estimators,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a set of elements to refine and coarsen.

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)
Template Parameters
ContainerContainer type for refinement and coarsen flags
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
Parameters
refineContainer of indices of elements to refine
coarsenContainer of indices of elements to coarsen
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [7/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( Func  flagFunction,
const Estimators &  estimators,
const DiscreteFuncs &  funcs,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh 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).
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
DiscreteFuncsType of the list of functions, must be random access container.
Parameters
flagFunctionFunction which takes an element index in range [0, Mesh()->ElementCount()), and returns an AdaptionType flag value indicating if the element should be refined, coarsened, or left 'as-is'
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
funcsList of functions defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [8/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( Func  flagFunction,
const Estimators &  estimators,
const DiscreteFunction< DIM, X, N, N > &  func,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh 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).
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
Parameters
flagFunctionFunction which takes an element index in range [0, Mesh()->ElementCount()), and returns an AdaptionType flag value indicating if the element should be refined, coarsened, or left 'as-is'
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
funcFunction defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ Adapt() [9/9]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::Adapt ( Func  flagFunction,
const Estimators &  estimators,
MeshSmoother  smoother = MeshSmoother{} 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh 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).
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
Parameters
flagFunctionFunction which takes an element index in range [0, Mesh()->ElementCount()), and returns an AdaptionType flag value indicating if the element should be refined, coarsened, or left 'as-is'
estimatorsError indicators for each element in mesh, necessary for checking convergence against expected
smootherSpecifies options for smoothing the mesh during refinemnt. Default options are as for a default constructed MeshSmoother object
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ AdaptConformingWithStrategy() [1/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptConformingWithStrategy ( const Estimators &  estimators,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor, with mesh smoothing options for conforming meshes (i.e., MeshSmoother::Conforming()).

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, estimators, MeshSmoother::Conforming())
FEMesh< DIM >::Modifications Adapt(const Container &flags, const Estimators &estimators, MeshSmoother smoother=MeshSmoother{})
Perform -adaptive refinement on the contained space and using a container of flags indicating the ref...
Definition: hpstrategies.hpp:502
static constexpr MeshSmoother Conforming(bool refineUnrefined=true)
Get a mesh smoother to create a conforming mesh with green refinement as necessary.
Definition: mesh.hpp:273
Exceptions
std::out_of_rangeMay be thrown if length of the estimators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::Conforming

◆ AdaptConformingWithStrategy() [2/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptConformingWithStrategy ( const Estimators &  estimators,
const DiscreteFuncs &  funcs,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified functions defined on the space with mesh smoothing options for conforming meshes (i.e., MeshSmoother::Conforming()).

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, estimators, funcs, MeshSmoother::Conforming())
Exceptions
std::out_of_rangeMay be thrown if length of the estimators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
DiscreteFuncsType of the list of functions, must be random access container.
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
funcsList of functions defined on the space to use for the refinement strategy
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::Conforming

◆ AdaptConformingWithStrategy() [3/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptConformingWithStrategy ( const Estimators &  estimators,
const DiscreteFunction< DIM, X, N, N > &  func,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified function defined on the space with mesh smoothing options for conforming meshes (i.e., MeshSmoother::Conforming()).

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, estimators, func, MeshSmoother::Conforming())
Exceptions
std::out_of_rangeMay be thrown if length of the estimators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
funcFunction defined on the space to use for the refinement strategy
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::Conforming

◆ AdaptWithStrategy() [1/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptWithStrategy ( const Estimators &  estimators,
const DiscreteFuncs &  funcs,
MeshSmoother  smoother,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified functions defined on the space and with the specified mesh smoothing options.

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, funcs, estimators, smoother)
Exceptions
std::out_of_rangeMay be thrown if length of the esitmators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
DiscreteFuncsType of the list of functions, must be random access container.
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
funcsList of functions defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt.
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)

◆ AdaptWithStrategy() [2/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptWithStrategy ( const Estimators &  estimators,
const DiscreteFunction< DIM, X, N, N > &  func,
MeshSmoother  smoother,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified function defined on the space and with the specified mesh smoothing options.

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, func, estimators, smoother)
Exceptions
std::out_of_rangeMay be thrown if length of the esitmators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
funcFunction defined on the space to use for the refinement strategy
smootherSpecifies options for smoothing the mesh during refinemnt.
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)

◆ AdaptWithStrategy() [3/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptWithStrategy ( const Estimators &  estimators,
MeshSmoother  smoother,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor and with the specified mesh smoothing options.

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, estimators, smoother)
Exceptions
std::out_of_rangeMay be thrown if length of the esitmators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
smootherSpecifies options for smoothing the mesh during refinemnt.
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)

◆ AdaptWithStrategyAndDefaultSmoothing() [1/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptWithStrategyAndDefaultSmoothing ( const Estimators &  estimators,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor, and with default mesh smoothing options (default constructed MeshSmoother).

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, estimators, MeshSmoother{})
Exceptions
std::out_of_rangeMay be thrown if length of the estimators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ AdaptWithStrategyAndDefaultSmoothing() [2/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptWithStrategyAndDefaultSmoothing ( const Estimators &  estimators,
const DiscreteFuncs &  funcs,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified functions defined on the space, and with default mesh smoothing options (default constructed MeshSmoother).

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, estimators, funcs, MeshSmoother{})
Exceptions
std::out_of_rangeMay be thrown if length of the estimators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
DiscreteFuncsType of the list of functions, must be random access container.
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
funcsList of functions defined on the space to use for the refinement strategy
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ AdaptWithStrategyAndDefaultSmoothing() [3/3]

FEMesh<DIM>::Modifications ptems::HPRefinementStrategy< DIM, double , N >::AdaptWithStrategyAndDefaultSmoothing ( const Estimators &  estimators,
const DiscreteFunction< DIM, X, N, N > &  func,
Args...  args 
)
inlineinherited

Perform \(hp\)-adaptive refinement on the contained space and mesh using a marking strategy functor object which takes error indicators as the first argument to the constructor with the specified function defined on the space, and with default mesh smoothing options (default constructed MeshSmoother).

Convenience method for call to

Adapt(MarkingStrategy{estimators, args...}, estimators, func, MeshSmoother{})
Exceptions
std::out_of_rangeMay be thrown if length of the estimators is not number of elements in mesh
Template Parameters
MarkingStrategyThe type of the marking strategy to use. The constructor of this type must take the error indicators as the first argument Must be explicitly specified
EstimatorsType of the list of error indicators, must be random access container of the same length as the number of elements on the mesh
ArgsType of extra arguments to pass to the MarkingStrategy constructor
Parameters
estimatorsError indicators for each element in mesh, passed to first argument of MarkingStrategy constructor.
funcFunction defined on the space to use for the refinement strategy
argsExtra arguments to pass to the MarkingStrategy constructor
Returns
Structure denoting changes that occurred to the mesh (but NOT the space polynomial degrees)
See also
MeshSmoother::MeshSmoother

◆ PostProcessMeshRefine()

virtual void ptems::HPRefinementStrategy< DIM, double , N >::PostProcessMeshRefine ( [[maybe_unused] ] const HPRefinementData< DIM, X, N > &  data,
[[maybe_unused] ] const typename FEMesh< DIM >::Modifications &  changes,
[[maybe_unused] ] const std::vector< AdaptationType > &  pRefinements 
)
inlineprotectedvirtualinherited

Method called by the adaptation routines after both \(h\)-refinement and \(p\)-refinement.

Can be used to for any post-processing.

Parameters
dataError indicators and functions defined on the mesh
changesModifications performed to the mesh ( \(h\)-refinements)
pRefinementsVector, with one entry per element, denoting the type of \(p\)-adaptation performed (refinement, coarsening, or none)

◆ PostProcessSpaceRefine()

virtual void ptems::HPRefinementStrategy< DIM, double , N >::PostProcessSpaceRefine ( [[maybe_unused] ] const HPRefinementData< DIM, X, N > &  data,
[[maybe_unused] ] const std::vector< AdaptationType > &  pRefinements 
)
inlineprotectedvirtualinherited

Method called by the adaptation routines after \(p\)-refinement, but before \(h\)-refinement.

Can be used to for any post-processing.

Parameters
dataError indicators and functions defined on the mesh
pRefinementsVector, with one entry per element, denoting the type of \(p\)-adaptation performed (refinement, coarsening, or none)

◆ ShouldPolydegAdapt() [1/2]

virtual bool ptems::HPRefinementStrategy< DIM, double , N >::ShouldPolydegAdapt ( const HPRefinementData< DIM, X, N > &  data,
std::size_t  element,
bool  refine 
)
protectedpure virtualinherited

Method called by the adaptation routines to decide if \(p\)-refinement should be performed for an element.

Parameters
dataError indicators and functions defiend on the mesh
elementIndex of the element
refinetrue if performing refinement of element, false if performing coarsening
Returns
true to perform \(p\)-refinement (or derefinement); false to perform \(h\)-refinement

◆ ShouldPolydegAdapt() [2/2]

template<std::size_t DIM, typename X = double, std::size_t N = 1>
virtual bool ptems::HPLegendreDecayStrategy< DIM, X, N >::ShouldPolydegAdapt ( const HPRefinementData< DIM, X, N > &  data,
std::size_t  element,
bool  refine 
)
inlineoverrideprotectedvirtual

Method called by the adaptation routines to decide if \(p\)-refinement should be performed for an element.

Parameters
dataError indicators and functions defiend on the mesh
elementIndex of the element
refinetrue if performing refinement of element, false if performing coarsening
Returns
true to perform \(p\)-refinement (or derefinement); false to perform \(h\)-refinement

◆ Space()

const PDiscreteFunctionSpace<DIM, double , N>& ptems::HPRefinementStrategy< DIM, double , N >::Space ( )
inlineinherited

Gets the finite element space the refinement strategy is on.


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