QUGaR 0.0.4
Loading...
Searching...
No Matches
qugar::impl Namespace Reference

Namespaces

namespace  funcs
 
namespace  tpms
 

Classes

class  AffineTransf
 Class for representing affine transformations. More...
 
class  BezierTP
 dim-dimensional tensor-product Bezier polynomial function. More...
 
class  DomainFunc
 Domain functions. More...
 
class  ImplReparamMesh
 Class for storing an implicit domain reparameterization using Lagrange cells. More...
 
class  MonomialsTP
 dim-dimensional tensor-product monomials function. More...
 
class  PolynomialTP
 Base class for tensor-product polynomial functions. More...
 
class  RefSystem
 A class representing a reference system in a given dimension. More...
 
struct  RootsIntervals
 Struct for storing and managing computed roots and intervals of an implicit function. More...
 
class  UnfittedImplDomain
 

Typedefs

template<int dim>
using ScalarFunc = DomainFunc<dim, 1>
 Alias for scalar functions.
 
template<int dim>
using ImplicitFunc = ScalarFunc<dim>
 Alias for implicit functions.
 

Enumerations

enum  FuncSign : std::int8_t { negative , positive , undetermined }
 

Functions

std::pair< Point< 3 >, Point< 3 > > create_reference_system_around_axis (Point< 3 > axis_z)
 
template<int dim, int range>
static bool is_bezier (const DomainFunc< dim, range > &func)
 Checks if a given DomainFunc object is of type BezierTP.
 
template<typename T >
void evaluate_Bernstein_value (const T &point, const int order, std::vector< T > &values)
 Evaluates the Bernstein polynomials of the given order.
 
template<typename T >
void evaluate_Bernstein (const T &point, const int order, int der, std::vector< T > &values)
 Evaluates the Bernstein polynomials of the given order (or its derivative).
 
template<int dim, int range>
std::shared_ptr< BezierTP< dim, range > > Bezier_product (const BezierTP< dim, range > &lhs, const BezierTP< dim, range > &rhs)
 Product of two Beziers.
 
template<int dim, int range, int dim2>
std::shared_ptr< BezierTP< dim, range > > Bezier_composition (const BezierTP< dim2, range > &lhs, const BezierTP< dim, dim2 > &rhs)
 Computes the composition of two Beziers as rhs(lhs)
 
template<int dim>
std::shared_ptr< const CutCellsQuad< dim > > create_quadrature (const UnfittedImplDomain< dim > &unf_domain, const std::vector< int > &cells, int n_pts_dir)
 
template<int dim>
std::shared_ptr< const CutUnfBoundsQuad< dim > > create_unfitted_bound_quadrature (const UnfittedImplDomain< dim > &unf_domain, const std::vector< int > &cells, int n_pts_dir)
 
template<int dim>
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > create_facets_quadrature (const UnfittedImplDomain< dim > &unf_domain, const std::vector< int > &cells, const std::vector< int > &facets, int n_pts_dir)
 
template<int dim, bool levelset>
std::shared_ptr< const ImplReparamMesh< levelset ? dim - 1 :dim, dim > > create_reparameterization (const UnfittedImplDomain< dim > &unf_domain, int n_pts_dir)
 
template<int dim, bool levelset>
std::shared_ptr< const ImplReparamMesh< levelset ? dim - 1 :dim, dim > > create_reparameterization (const UnfittedImplDomain< dim > &unf_domain, const std::vector< int > &cells, int n_pts_dir)
 
template<int dim, bool S = false>
std::shared_ptr< ImplReparamMesh< S ? dim - 1 :dim, dim > > reparam_Bezier (const BezierTP< dim, 1 > &bzr, const BoundBox< dim > &domain, int order)
 Reparameterizes a Bezier implicit function in the unit hypercube domain.
 
template<int dim, bool S = false>
void reparam_Bezier (const BezierTP< dim, 1 > &bzr, const BoundBox< dim > &domain, ImplReparamMesh< S ? dim - 1 :dim, dim > &reparam)
 Reparameterizes a Bezier implicit function in the unit hypercube domain.
 
template<int dim, bool S = false>
std::shared_ptr< ImplReparamMesh< S ? dim - 1 :dim, dim > > reparam_general (const ImplicitFunc< dim > &func, const BoundBox< dim > &domain, int order)
 Reparameterizes the domain defined by a general implicit function.
 
template<int dim, bool S = false>
void reparam_general (const ImplicitFunc< dim > &func, const BoundBox< dim > &domain, ImplReparamMesh< S ? dim - 1 :dim, dim > &reparam)
 Reparameterizes the domain defined by a general implicit function.
 
template<int dim>
std::shared_ptr< ImplReparamMesh< dim - 1, dim > > reparam_general_facet (const ImplicitFunc< dim > &func, const BoundBox< dim > &domain, int facet_id, int order)
 Reparameterizes a face of domain defined by an implicit function. It reparameterizes one of the 2*dim faces of the domain.
 
template<int dim>
void reparam_general_facet (const ImplicitFunc< dim > &func, const BoundBox< dim > &domain, int facet_id, ImplReparamMesh< dim - 1, dim > &reparam)
 Reparameterizes a face of domain defined by an implicit function. It reparameterizes one of the 2*dim faces of the domain.
 
template<int dim>
bool on_levelset (const ImplicitFunc< dim > &phi, const Point< dim > &point, Tolerance tol=Tolerance())
 Checks if a point belongs to the levelset of an implicit function phi.
 
template<int dim>
int get_facet_constant_dir (int local_facet_id)
 Gets the constant direction of the local facet.
 
template<int dim>
int get_facet_side (int local_facet_id)
 Gets the side of the facet. Either 0 or 1.
 
template<int dim>
int get_local_facet_id (int const_dir, int side)
 Get the local facet ID for a given const direction and side.
 
template<int dim>
Vector< int, dim - 1 > get_edge_constant_dirs (int local_edge_id)
 Gets the constant directions of the local edge.
 
template<int dim>
Vector< int, dim - 1 > get_edge_sides (int local_edge_id)
 Gets the sides of the edge. Either 0 or 1 along each constant direction.
 
void evaluate_Lagrange_basis_1D (real point, int order, bool chebyshev, std::vector< real > &values)
 Evaluates the Lagrange basis polynomials in 1D at a given point.
 
void evaluate_Lagrange_basis_der_1D (real point, int order, bool chebyshev, std::vector< real > &values)
 Evaluates the first derivative of the Lagrange basis polynomial in 1D.
 
template<int dim>
void evaluate_Lagrange_basis (const Point< dim > &point, const TensorSizeTP< dim > &order, bool chebyshev, std::vector< real > &basis)
 Evaluates the tensor-product Lagrange basis functions at a given point.
 
template<int dim>
void evaluate_Lagrange_derivative (const Point< dim > &point, const TensorSizeTP< dim > &order, bool chebyshev, Vector< std::vector< real >, dim > &basis_ders)
 Evaluates the derivative of the Lagrange basis functions at a given point, along all directions.
 

Typedef Documentation

◆ ImplicitFunc

template<int dim>
using qugar::impl::ImplicitFunc = ScalarFunc<dim>

Alias for implicit functions.

◆ ScalarFunc

template<int dim>
using qugar::impl::ScalarFunc = DomainFunc<dim, 1>

Alias for scalar functions.

Enumeration Type Documentation

◆ FuncSign

enum qugar::impl::FuncSign : std::int8_t
Enumerator
negative 
positive 
undetermined 

Function Documentation

◆ Bezier_composition()

template<int dim, int range, int dim2>
std::shared_ptr< BezierTP< dim, range > > qugar::impl::Bezier_composition ( const BezierTP< dim2, range > & lhs,
const BezierTP< dim, dim2 > & rhs )
nodiscard

Computes the composition of two Beziers as rhs(lhs)

Template Parameters
dimDimension of the resultant Bezier.
rangeRange of the resultant Bezier.
dim2Common dimension between Beziers.
Parameters
lhsBezier to be composed.
rhsComposing Bezier.
Returns
Resultant composition.

◆ Bezier_product()

template<int dim, int range>
std::shared_ptr< BezierTP< dim, range > > qugar::impl::Bezier_product ( const BezierTP< dim, range > & lhs,
const BezierTP< dim, range > & rhs )
nodiscard

Product of two Beziers.

Parameters
lhsFirst Bezier to multiply.
rhsSecond Bezier to multiply.
Returns
Product of Beziers.

◆ create_facets_quadrature()

template<int dim>
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > qugar::impl::create_facets_quadrature ( const UnfittedImplDomain< dim > & unf_domain,
const std::vector< int > & cells,
const std::vector< int > & facets,
int n_pts_dir )

◆ create_quadrature()

template<int dim>
std::shared_ptr< const CutCellsQuad< dim > > qugar::impl::create_quadrature ( const UnfittedImplDomain< dim > & unf_domain,
const std::vector< int > & cells,
int n_pts_dir )
Here is the caller graph for this function:

◆ create_reference_system_around_axis()

std::pair< Point< 3 >, Point< 3 > > qugar::impl::create_reference_system_around_axis ( Point< 3 > axis_z)

◆ create_reparameterization() [1/2]

template<int dim, bool levelset>
std::shared_ptr< const ImplReparamMesh< levelset ? dim - 1 :dim, dim > > qugar::impl::create_reparameterization ( const UnfittedImplDomain< dim > & unf_domain,
const std::vector< int > & cells,
int n_pts_dir )

◆ create_reparameterization() [2/2]

template<int dim, bool levelset>
std::shared_ptr< const ImplReparamMesh< levelset ? dim - 1 :dim, dim > > qugar::impl::create_reparameterization ( const UnfittedImplDomain< dim > & unf_domain,
int n_pts_dir )

◆ create_unfitted_bound_quadrature()

template<int dim>
std::shared_ptr< const CutUnfBoundsQuad< dim > > qugar::impl::create_unfitted_bound_quadrature ( const UnfittedImplDomain< dim > & unf_domain,
const std::vector< int > & cells,
int n_pts_dir )
Here is the caller graph for this function:

◆ evaluate_Bernstein()

template<typename T >
void qugar::impl::evaluate_Bernstein ( const T & point,
const int order,
int der,
std::vector< T > & values )

Evaluates the Bernstein polynomials of the given order (or its derivative).

Template Parameters
TType of the input coordinate.
VType of the output vector of basis values.
Parameters
pointEvaluation point.
orderOrder of the polynomial.
derOrder of the derivative to be computed. If 0, the value itself is computed.
valuesComputed basis values.

◆ evaluate_Bernstein_value()

template<typename T >
void qugar::impl::evaluate_Bernstein_value ( const T & point,
const int order,
std::vector< T > & values )

Evaluates the Bernstein polynomials of the given order.

Template Parameters
TType of the input coordinate.
Parameters
pointEvaluation point.
orderOrder of the polynomial.
valuesComputed basis values.

◆ evaluate_Lagrange_basis()

template<int dim>
void qugar::impl::evaluate_Lagrange_basis ( const Point< dim > & point,
const TensorSizeTP< dim > & order,
bool chebyshev,
std::vector< real > & basis )

Evaluates the tensor-product Lagrange basis functions at a given point.

Parameters
pointThe point at which to evaluate the Lagrange basis functions.
orderThe order (degree + 1) of the tensor product along each direction.
chebyshevA boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).
basisA reference to a vector where the computed derivative values will be stored. This vector will be resized to the order of the Lagrange basis polynomials.

◆ evaluate_Lagrange_basis_1D()

void qugar::impl::evaluate_Lagrange_basis_1D ( real point,
int order,
bool chebyshev,
std::vector< real > & values )

Evaluates the Lagrange basis polynomials in 1D at a given point.

This function computes the values of the Lagrange basis polynomials of a specified order at a given point. The basis can be evaluated using either Chebyshev nodes (of 2nd kind) or standard equidistant nodes.

Parameters
pointThe point at which to evaluate the Lagrange basis polynomials.
orderThe order of the Lagrange basis polynomials (degree + 1).
chebyshevA boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).
valuesA reference to a vector where the computed values of the Lagrange basis polynomials will be stored. This vector will be resized to the order of the Lagrange basis polynomials.

◆ evaluate_Lagrange_basis_der_1D()

void qugar::impl::evaluate_Lagrange_basis_der_1D ( real point,
int order,
bool chebyshev,
std::vector< real > & values )

Evaluates the first derivative of the Lagrange basis polynomial in 1D.

This function computes the values of the first derivative of the Lagrange basis polynomial at a given point for a specified order. The basis can be either Chebyshev (of 2nd kind) or equidistant nodes.

Parameters
pointThe point at which to evaluate the derivative.
orderThe order of the Lagrange basis polynomial (degree + 1).
chebyshevA boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).
valuesA reference to a vector where the computed derivative values will be stored. This vector will be resized to the order of the Lagrange basis polynomials.

◆ evaluate_Lagrange_derivative()

template<int dim>
void qugar::impl::evaluate_Lagrange_derivative ( const Point< dim > & point,
const TensorSizeTP< dim > & order,
bool chebyshev,
Vector< std::vector< real >, dim > & basis_ders )

Evaluates the derivative of the Lagrange basis functions at a given point, along all directions.

Parameters
pointThe point at which to evaluate the derivative of the Lagrange basis functions.
orderThe order (degree + 1) of the tensor product along each direction.
chebyshevA boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).
basis_dersA vector to store the computed derivatives of the Lagrange basis functions along all directions. The vectors along each direction will be resized to the order of the Lagrange basis polynomials.

◆ get_edge_constant_dirs()

template<int dim>
Vector< int, dim - 1 > qugar::impl::get_edge_constant_dirs ( int local_edge_id)
nodiscard

Gets the constant directions of the local edge.

The edges are assumed to follow a lexicographical numbering.

Parameters
local_edge_idId of the edge referred to a cell.
Returns
Constant directions in the range [0,dim[.

◆ get_edge_sides()

template<int dim>
Vector< int, dim - 1 > qugar::impl::get_edge_sides ( int local_edge_id)
nodiscard

Gets the sides of the edge. Either 0 or 1 along each constant direction.

Along a constant direction, it returns if the local edge corresponds to the side of the cell's bounding box with minimum coordinate along the constant direction of the facet (0), or the maximum (1).

Parameters
local_edge_idId of the edge referred to a cell.
Returns
Sides of the edge.

◆ get_facet_constant_dir()

template<int dim>
int qugar::impl::get_facet_constant_dir ( int local_facet_id)
nodiscard

Gets the constant direction of the local facet.

The constant direction is computed as the floor of division by 2 of local_facet_id.

Parameters
local_facet_idId of the facet referred to a cell. It must be a value in the range [0, 2*dim[.
Returns
Constant direction in the range [0,dim[.

◆ get_facet_side()

template<int dim>
int qugar::impl::get_facet_side ( int local_facet_id)
nodiscard

Gets the side of the facet. Either 0 or 1.

I.e., it returns if the local facet corresponds to the side of the cell's bounding box with minimum coordinate along the constant direction of the facet (0), or the maximum (1).

Side 0 correspond to even values of local_facet_id, and side 1 to odd values.

Parameters
local_facet_idId of the facet referred to a cell. It must be a value in the range [0, 2*dim[.
Returns
Side of the facet. Either 0 or 1.

◆ get_local_facet_id()

template<int dim>
int qugar::impl::get_local_facet_id ( int const_dir,
int side )
nodiscard

Get the local facet ID for a given const direction and side.

This function returns the local facet ID based on the specified constant direction and side. It is used in the context of a multi-dimensional space where facets (or faces) of a geometric entity are identified.

Template Parameters
dimThe dimension of the space.
Parameters
const_dirThe constant direction for which the local facet ID is to be determined.
sideThe side (0 or 1) in the specified direction.
Returns
The local facet ID corresponding to the given direction and side.

◆ is_bezier()

template<int dim, int range>
static bool qugar::impl::is_bezier ( const DomainFunc< dim, range > & func)
static

Checks if a given DomainFunc object is of type BezierTP.

This function uses dynamic_cast to determine if the provided DomainFunc object can be cast to a BezierTP object. If the cast is successful, the function returns true, indicating that the object is of type BezierTP. Otherwise, it returns false.

Template Parameters
dimThe dimension of the BezierTP object.
rangeThe range of the BezierTP object.
Parameters
funcThe DomainFunc object to be checked.
Returns
true if the object is of type BezierTP, false otherwise.

◆ on_levelset()

template<int dim>
bool qugar::impl::on_levelset ( const ImplicitFunc< dim > & phi,
const Point< dim > & point,
Tolerance tol = Tolerance() )
nodiscard

Checks if a point belongs to the levelset of an implicit function phi.

Template Parameters
dimDimension of the function and point.
Parameters
phiImplicit function to be tested.
pointPoint to be tested.
tolTolerance to be used for checking if the value of phi is (close to) zero. If not provided, default tolerance is used.
Returns
Whether the point belong to the levelset of phi.

◆ reparam_Bezier() [1/2]

template<int dim, bool S = false>
void qugar::impl::reparam_Bezier ( const BezierTP< dim, 1 > & bzr,
const BoundBox< dim > & domain,
ImplReparamMesh< S ? dim - 1 :dim, dim > & reparam )

Reparameterizes a Bezier implicit function in the unit hypercube domain.

Note
The generated reparameterization has a wirebasket associated, but coincident were not merged.
Template Parameters
dimParametric dimension of the function.
SFlag indicating if the reparameterization must be performed only for the levelset surface (true), i.e., the manifold where the Bezier function is equal to 0, or the volume (false), i.e., the subregion where the Bezier function is negative.
Parameters
bzrBezier implicit function to reparameterize.
domainDomain to which the implicit function refers to (even if Beziers are defined in the unit domain).
reparamReparameterization container to which new generated cells are appended to.

◆ reparam_Bezier() [2/2]

template<int dim, bool S = false>
std::shared_ptr< ImplReparamMesh< S ? dim - 1 :dim, dim > > qugar::impl::reparam_Bezier ( const BezierTP< dim, 1 > & bzr,
const BoundBox< dim > & domain,
int order )

Reparameterizes a Bezier implicit function in the unit hypercube domain.

Note
The generated reparameterization has a wirebasket associated, but coincident were not merged.
Template Parameters
dimParametric dimension of the function.
SFlag indicating if the reparameterization must be performed only for the levelset surface (true), i.e., the manifold where the Bezier function is equal to 0, or the volume (false), i.e., the subregion where the Bezier function is negative.
Parameters
bzrBezier implicit function to reparameterize.
domainDomain to which the implicit function refers to (even if Beziers are defined in the unit domain).
orderOrder of the reparameterization (number of points per direction in each reparameterization cell).
Returns
Generated reparameterization.

◆ reparam_general() [1/2]

template<int dim, bool S = false>
void qugar::impl::reparam_general ( const ImplicitFunc< dim > & func,
const BoundBox< dim > & domain,
ImplReparamMesh< S ? dim - 1 :dim, dim > & reparam )

Reparameterizes the domain defined by a general implicit function.

Note
The generated reparameterization has a wirebasket associated, but coincident were not merged.
Template Parameters
dimParametric dimension of the function.
SFlag indicating if the reparameterization must be performed only for the levelset surface (true), i.e., the manifold where the either of the two Bezier functions are equal to 0, or the subregion (false) between those surfaces, where both Bezier functions are negative.
Parameters
funcFunction to reparameterize.
domainDomain to reparameterize.
reparamReparameterization container to which new generated cells are appended to.

◆ reparam_general() [2/2]

template<int dim, bool S = false>
std::shared_ptr< ImplReparamMesh< S ? dim - 1 :dim, dim > > qugar::impl::reparam_general ( const ImplicitFunc< dim > & func,
const BoundBox< dim > & domain,
int order )

Reparameterizes the domain defined by a general implicit function.

Note
The generated reparameterization has a wirebasket associated, but coincident were not merged.
Template Parameters
dimParametric dimension of the function.
SFlag indicating if the reparameterization must be performed only for the levelset surface (true), i.e., the manifold where the either of the two Bezier functions are equal to 0, or the subregion (false) between those surfaces, where both Bezier functions are negative.
Parameters
funcFunction to reparameterize.
domainDomain to reparameterize.
orderOrder of the reparameterization (number of points per direction in each reparameterization cell).
Returns
Generated reparameterization.

◆ reparam_general_facet() [1/2]

template<int dim>
void qugar::impl::reparam_general_facet ( const ImplicitFunc< dim > & func,
const BoundBox< dim > & domain,
int facet_id,
ImplReparamMesh< dim - 1, dim > & reparam )

Reparameterizes a face of domain defined by an implicit function. It reparameterizes one of the 2*dim faces of the domain.

Note
The generated reparameterization has a wirebasket associated, but coincident were not merged.
Template Parameters
dimParametric dimension of the function.
Parameters
funcFunction to reparameterize.
domainDomain to reparameterize.
facet_idId of the face to reparameterize. It must be a value in the range [0, 2*dim[.
reparamReparameterization container to which new generated cells are appended to.

◆ reparam_general_facet() [2/2]

template<int dim>
std::shared_ptr< ImplReparamMesh< dim - 1, dim > > qugar::impl::reparam_general_facet ( const ImplicitFunc< dim > & func,
const BoundBox< dim > & domain,
int facet_id,
int order )

Reparameterizes a face of domain defined by an implicit function. It reparameterizes one of the 2*dim faces of the domain.

Note
The generated reparameterization has a wirebasket associated, but coincident were not merged.
Template Parameters
dimParametric dimension of the function.
Parameters
funcFunction to reparameterize.
domainDomain to reparameterize.
facet_idId of the face to reparameterize. It must be a value in the range [0, 2*dim[.
orderOrder of the reparameterization (number of points per direction in each reparameterization cell).
Returns
Generated reparameterization.