QUGaR 0.0.9
Loading...
Searching...
No Matches
qugar Namespace Reference

QUGaR's main namespace. More...

Namespaces

namespace  impl
 
namespace  numbers
 Numbers' namespace.
 

Classes

class  BoundBox
 Class representing a dim-dimensional Cartesian product bounding box. More...
 
class  CartGridTP
 Class representing a dim-dimensional Cartesian tensor-product grid. More...
 
struct  CutCellsQuad
 
struct  CutIsoBoundsQuad
 
struct  CutUnfBoundsQuad
 
class  Gauss
 Helper class for computing the abscissae and weights of the Gauss-Legendre quadrature rule referred to the [0,1]. The maximum number of abscissae/weights is prescribed by n_max. The abscissae / weights are precomputed with enough precision (34 significant digits) for double or quadruple precision applications. More...
 
class  Quadrature
 Class for storing dim-dimensional quadratures (non-tensor product). More...
 
class  ReparamMesh
 Class for storing an implicit domain reparameterization using Lagrange cells. More...
 
class  SubCartGridTP
 Subgrid of a Cartesian grid TP. It is a subset of the cells of a given grid. More...
 
class  TanhSinh
 Helper class for computing the abscissae and weights of the tanh-sinh quadrature rule referred to the [0,1]. The maximum number of abscissae/weights is prescribed by n_max. The abscissae / weights are precomputed with enough precision (34 significant digits) for double or quadruple precision applications. More...
 
class  TensorIndexRangeTP
 Class representing a dim-dimensional range defined by lower and upper tensor bounds. More...
 
class  TensorIndexTP
 Class representing a dim-dimensional tensor-product indices. More...
 
class  TensorSizeTP
 Class representing a dim-dimensional tensor-product sizes container. More...
 
class  Tolerance
 Class for tolerance related computations. More...
 
class  UnfittedBinarySpacePart
 
class  UnfittedDomain
 

Concepts

concept  Is2Dor3D
 Checks if the given dimension is 2 or 3.
 
concept  Is1D
 Checks if the given dimension is 1.
 

Typedefs

using Interval = BoundBox<1>
 Alias representing an interval as a 1-dimensional bounding box.
 
template<int dim, typename T = real>
using Point = Vector<T, dim>
 Class representing a dim-dimensional Point.
 
using real = double
 
using index = std::ptrdiff_t
 
template<typename T , int dim>
using Vector = ::algoim::uvector<T, dim>
 Class representing a vector.
 

Enumerations

enum  ImmersedStatus : std::uint8_t { cut , full , empty }
 
enum class  ImmersedFacetStatus : std::uint8_t {
  cut , full , empty , cut_unf_bdry ,
  cut_ext_bdry , cut_unf_bdry_ext_bdry , full_unf_bdry , full_ext_bdry ,
  unf_bdry , ext_bdry , unf_bdry_ext_bdry
}
 
enum class  ImmersedCellStatus : std::uint8_t { cut , full , empty , unknown }
 

Functions

template<int dim>
std::shared_ptr< const CutCellsQuad< dim > > create_quadrature (const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, int n_pts_dir)
 Creates quadrature for cells.
 
template<int dim>
std::shared_ptr< const CutUnfBoundsQuad< dim > > create_unfitted_bound_quadrature (const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, int n_pts_dir, bool include_facet_unf_bdry, bool exclude_ext_bdry)
 Creates a quadrature for the unfitted boundary.
 
template<int dim>
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > create_facets_quadrature_interior_integral (const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, const std::vector< int > &facets, int n_pts_dir)
 Creates quadrature for interior integrals.
 
template<int dim>
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > create_facets_quadrature_exterior_integral (const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, const std::vector< int > &facets, int n_pts_dir)
 Creates quadrature for exterior integrals.
 
template<int dim>
void find_coincident_points (const std::vector< Point< dim > > &points, const Tolerance &tol, std::vector< std::size_t > &points_map)
 Finds coincident points in a given set of points.
 
template<int dim>
void make_points_unique (std::vector< Point< dim > > &points, const Tolerance &tol, std::vector< std::size_t > &old_to_new)
 Makes a given vector of points unique by merging coincident points.
 
template<int dim, bool levelset>
std::shared_ptr< const ReparamMesh< levelset ? dim - 1 :dim, dim > > create_reparameterization (const UnfittedDomain< dim > &unf_domain, int n_pts_dir, bool merge_points)
 
template<int dim, bool levelset>
std::shared_ptr< const ReparamMesh< levelset ? dim - 1 :dim, dim > > create_reparameterization (const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, int n_pts_dir, bool merge_points)
 
template<class T , class V >
constexpr T narrow_cast (V &&val) noexcept
 
template<class T , std::size_t N>
constexpr T & at (T(&arr)[N], const index ind)
 
template<class Cont >
constexpr auto at (Cont &cont, const index ind) -> decltype(cont[cont.size()])
 
template<class T >
constexpr T at (const std::initializer_list< T > cont, const index ind)
 
template<class T , std::size_t extent = std::dynamic_extent>
constexpr auto at (std::span< T, extent > spn, const index ind) -> decltype(spn[spn.size()])
 
template<typename T , int dim>
Vector< T, dim > permute_vector_directions (const Vector< T, dim > &vec)
 

Variables

constexpr int QUGAR_VERSION_MAJOR = { 0 }
 
constexpr int QUGAR_VERSION_MINOR = { 0 }
 
constexpr int QUGAR_VERSION_PATCH = { 9 }
 
constexpr int QUGAR_VERSION_TWEAK = { }
 
constexpr std::string_view QUGAR_VERSION_STRING { "0.0.9" }
 
constexpr std::string_view QUGAR_VERSION_GIT { "91692c2" }
 

Detailed Description

QUGaR's main namespace.

Typedef Documentation

◆ index

using qugar::index = std::ptrdiff_t

◆ Interval

Alias representing an interval as a 1-dimensional bounding box.

◆ Point

template<int dim, typename T = real>
using qugar::Point = Vector<T, dim>

Class representing a dim-dimensional Point.

Template Parameters
dimDimension of point.

◆ real

using qugar::real = double

◆ Vector

template<typename T , int dim>
using qugar::Vector = ::algoim::uvector<T, dim>

Class representing a vector.

Template Parameters
dimDimension of point.
TVector type.

Enumeration Type Documentation

◆ ImmersedCellStatus

enum class qugar::ImmersedCellStatus : std::uint8_t
strong
Enumerator
cut 
full 
empty 
unknown 

◆ ImmersedFacetStatus

enum class qugar::ImmersedFacetStatus : std::uint8_t
strong
Enumerator
cut 
full 
empty 
cut_unf_bdry 
cut_ext_bdry 
cut_unf_bdry_ext_bdry 
full_unf_bdry 
full_ext_bdry 
unf_bdry 
ext_bdry 
unf_bdry_ext_bdry 

◆ ImmersedStatus

enum qugar::ImmersedStatus : std::uint8_t
Enumerator
cut 
full 
empty 

Function Documentation

◆ at() [1/4]

template<class T >
T qugar::at ( const std::initializer_list< T > cont,
const index ind )
constexpr
Here is the call graph for this function:

◆ at() [2/4]

template<class Cont >
auto qugar::at ( Cont & cont,
const index ind ) -> decltype(cont[cont.size()])
constexpr
Here is the call graph for this function:

◆ at() [3/4]

template<class T , std::size_t extent = std::dynamic_extent>
auto qugar::at ( std::span< T, extent > spn,
const index ind ) -> decltype(spn[spn.size()])
constexpr
Here is the call graph for this function:

◆ at() [4/4]

template<class T , std::size_t N>
T & qugar::at ( T(&) arr[N],
const index ind )
constexpr
Here is the call graph for this function:

◆ create_facets_quadrature_exterior_integral()

template<int dim>
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > qugar::create_facets_quadrature_exterior_integral ( const UnfittedDomain< dim > & unf_domain,
const std::vector< std::int64_t > & cells,
const std::vector< int > & facets,
int n_pts_dir )

Creates quadrature for exterior integrals.

It creates quadrature for the active part of the facets that belong to the external boundary of the domain. The external boundary of the domain may be the external boundary of the grid or the unfitted boundary.

Template Parameters
dimDimension of the domain.
Parameters
unf_domainUnfitted domain.
cellsCells for which the quadrature is created.
facetsLocal facets for which the quadrature is created (this vector must be of the same size as cells).
n_pts_dirNumber of points in each direction for generated custom quadratures.
Returns
Generated quadratures.

◆ create_facets_quadrature_interior_integral()

template<int dim>
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > qugar::create_facets_quadrature_interior_integral ( const UnfittedDomain< dim > & unf_domain,
const std::vector< std::int64_t > & cells,
const std::vector< int > & facets,
int n_pts_dir )

Creates quadrature for interior integrals.

It creates quadratures for the internal part of interior facets. I.e., fully internal facets or the cut (interior) part of cut facets. Unfitted boundaries laying on the facet are not considered.

Warning
The provided facets are not checked to be interior facets. It is the caller's responsibility to provide the correct facets.
Template Parameters
dimDimension of the domain.
Parameters
unf_domainUnfitted domain.
cellsCells for which the quadrature is created.
facetsLocal facets for which the quadrature is created (this vector must be of the same size as cells).
n_pts_dirNumber of points in each direction for generated custom quadratures.
Returns
Generated quadratures.

◆ create_quadrature()

template<int dim>
std::shared_ptr< const CutCellsQuad< dim > > qugar::create_quadrature ( const UnfittedDomain< dim > & unf_domain,
const std::vector< std::int64_t > & cells,
int n_pts_dir )

Creates quadrature for cells.

For cut cells, the generated quadrature correspond to the interior part of the cell. For full cells, the quadrature is the standard one, and for empty cells, no quadrature is created (the generated quadrature is empty except for the number of points that is set to 0).

Template Parameters
dimDimension of the domain.
Parameters
unf_domainUnitted domain.
cellsCells for which the quadrature is created.
n_pts_dirNumber of points in each direction for generated custom quadratures.
Returns
Generated quadrature.
Here is the caller graph for this function:

◆ create_reparameterization() [1/2]

template<int dim, bool levelset>
std::shared_ptr< const ReparamMesh< levelset ? dim - 1 :dim, dim > > qugar::create_reparameterization ( const UnfittedDomain< dim > & unf_domain,
const std::vector< std::int64_t > & cells,
int n_pts_dir,
bool merge_points )

◆ create_reparameterization() [2/2]

template<int dim, bool levelset>
std::shared_ptr< const ReparamMesh< levelset ? dim - 1 :dim, dim > > qugar::create_reparameterization ( const UnfittedDomain< dim > & unf_domain,
int n_pts_dir,
bool merge_points )

◆ create_unfitted_bound_quadrature()

template<int dim>
std::shared_ptr< const CutUnfBoundsQuad< dim > > qugar::create_unfitted_bound_quadrature ( const UnfittedDomain< dim > & unf_domain,
const std::vector< std::int64_t > & cells,
int n_pts_dir,
bool include_facet_unf_bdry,
bool exclude_ext_bdry )

Creates a quadrature for the unfitted boundary.

Template Parameters
dimDimension of the domain.
Parameters
unf_domainUnitted domain.
cellsCells for which the quadrature is created.
n_pts_dirNumber of points in each direction for generated custom quadratures.
include_facet_unf_bdryIf true, the quadrature includes the parts of the unfitted boundary that belong to the cells' facets.
exclude_ext_bdryIf the previous parameter is true, and this is one is also true, the parts of the unfitted boundary that belong to the external facets are not included.
Returns
Generated quadrature.
Here is the caller graph for this function:

◆ find_coincident_points()

template<int dim>
void qugar::find_coincident_points ( const std::vector< Point< dim > > & points,
const Tolerance & tol,
std::vector< std::size_t > & points_map )

Finds coincident points in a given set of points.

This function identifies pairs of points in the provided vector that are coincident within a specified tolerance. It returns a vector that maps the (indices of the) points in the input vector, to the unique ones. If more than two points are coincident, all those points will be associated to a single one.

Template Parameters
dimThe dimension of the points.
Parameters
pointsA vector of points to be checked for coincidences.
tolThe tolerance within which points are considered coincident.
points_mapA vector mapping original point indices to unique ones.

◆ make_points_unique()

template<int dim>
void qugar::make_points_unique ( std::vector< Point< dim > > & points,
const Tolerance & tol,
std::vector< std::size_t > & old_to_new )

Makes a given vector of points unique by merging coincident points.

This function identifies pairs of points in the provided vector that are coincident within a specified tolerance, and merges them into a single one. It returns a vector that map the (indices of the) original points, to the new (unique) ones. If more than two points are coincident, all those points will be associated to a single one.

Template Parameters
dimThe dimension of the points.
Parameters
pointsVector of points to make unique.
tolThe tolerance within which points are considered coincident.
old_to_newResult vector mapping original point indices to the new (unique) ones.

◆ narrow_cast()

template<class T , class V >
T qugar::narrow_cast ( V && val)
constexprnoexcept
Here is the caller graph for this function:

◆ permute_vector_directions()

template<typename T , int dim>
Vector< T, dim > qugar::permute_vector_directions ( const Vector< T, dim > & vec)

Permutes the directions of a given vector by reversing the order of its elements.

Note
This function is useful for transforming points for some Algoim's algorithms that use counter-lexicographical ordering (as coefficients of Bezier polynomials).
Template Parameters
TThe type of the elements in the vector.
dimThe dimension of the vector.
Parameters
vecThe input vector to be permuted.
Returns
A new vector with its directions permuted.

Variable Documentation

◆ QUGAR_VERSION_GIT

std::string_view qugar::QUGAR_VERSION_GIT { "91692c2" }
constexpr

◆ QUGAR_VERSION_MAJOR

int qugar::QUGAR_VERSION_MAJOR = { 0 }
constexpr

◆ QUGAR_VERSION_MINOR

int qugar::QUGAR_VERSION_MINOR = { 0 }
constexpr

◆ QUGAR_VERSION_PATCH

int qugar::QUGAR_VERSION_PATCH = { 9 }
constexpr

◆ QUGAR_VERSION_STRING

std::string_view qugar::QUGAR_VERSION_STRING { "0.0.9" }
constexpr

◆ QUGAR_VERSION_TWEAK

int qugar::QUGAR_VERSION_TWEAK = { }
constexpr