C++

This is experimental documentation for the C++ API. Here below you will find an unstructured listing of all the available functions and classes in the C++ API. Feel free to visit the contributors section is you want to help!

The full Doxygen generated documentation can be visited here.

namespace qugar

QUGaR’s main namespace.

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.

Template Parameters:

dim – Dimension of point.

using real = double
using index = std::ptrdiff_t
template<typename T, int dim>
using Vector = ::algoim::uvector<T, dim>

Class representing a vector.

Template Parameters:
  • dim – Dimension of point.

  • T – Vector type.

Enums

enum ImmersedStatus

Values:

enumerator cut
enumerator full
enumerator empty
enum class ImmersedFacetStatus : std::uint8_t

Values:

enumerator cut
enumerator full
enumerator empty
enumerator cut_unf_bdry
enumerator cut_ext_bdry
enumerator cut_unf_bdry_ext_bdry
enumerator full_unf_bdry
enumerator full_ext_bdry
enumerator unf_bdry
enumerator ext_bdry
enumerator unf_bdry_ext_bdry

Functions

template<int dim>
std::shared_ptr<const CutCellsQuad<dim>> create_quadrature(const UnfittedDomain<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 UnfittedDomain<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 UnfittedDomain<dim> &unf_domain, const std::vector<int> &cells, const std::vector<int> &facets, int n_pts_dir)
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.

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:

dim – The dimension of the points.

Parameters:
  • points – A vector of points to be checked for coincidences.

  • tol – The tolerance within which points are considered coincident.

  • points_map – A vector mapping original point indices to unique ones.

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.

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:

dim – The dimension of the points.

Parameters:
  • points – Vector of points to make unique.

  • tol – The tolerance within which points are considered coincident.

  • old_to_new – Result vector mapping original point indices to the new (unique) ones.

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)
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<int> &cells, int n_pts_dir)
template<class T, class V>
T narrow_cast(V &&val) noexcept
template<class T, std::size_t N>
T &at(T (&arr)[N], const index ind)
template<class Cont>
decltype(cont[cont.size()]) at(Cont &cont, const index ind)
template<class T>
T at(const std::initializer_list<T> cont, const index ind)
template<class T, std::size_t extent = std::dynamic_extent>
decltype(spn[spn.size()]) at(std::span<T, extent> spn, const index ind)
template<typename T, int dim>
Vector<T, dim> 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:
  • T – The type of the elements in the vector.

  • dim – The dimension of the vector.

Parameters:

vec – The input vector to be permuted.

Returns:

A new vector with its directions permuted.

Variables

int QUGAR_VERSION_MAJOR = {0}
int QUGAR_VERSION_MINOR = {0}
int QUGAR_VERSION_PATCH = {4}
int QUGAR_VERSION_TWEAK = {}
std::string_view QUGAR_VERSION_STRING = {"0.0.4"}
std::string_view QUGAR_VERSION_GIT = {"c9cb306"}
template<int dim>
class BoundBox
#include <bbox.hpp>

Class representing a dim-dimensional Cartesian product bounding box.

Template Parameters:

dim – Dimension of the domain.

Constructors

BoundBox()

Construct a new default BoundBox. Initializes the box to the unit_ domain [0, 1].

BoundBox(const std::array<real, static_cast<std::size_t>(dim)> &min, const std::array<real, static_cast<std::size_t>(dim)> &max)

Construct a new BoundBox object from the max and min coordinates.

Parameters:
  • min – Minimum coordinates of the box.

  • max – Maximum coordinates of the box.

BoundBox(const real min, const real max)

Construct a new BoundBox object from the max and min coordinates.

Parameters:
  • min – Minimum coordinates of the box in all directions.

  • max – Maximum coordinates of the box in all directions.

BoundBox(const Point<dim> &min, const Point<dim> &max)

Construct a new BoundBox object from the max and min coordinates.

Parameters:
  • min – Minimum coordinates of the box.

  • max – Maximum coordinates of the box.

explicit BoundBox(const ::algoim::HyperRectangle<real, dim> &rectangle)

Constructs a BoundBox from a given Algoim’s hyperrectangle.

This constructor initializes a BoundBox object using the provided hyperrectangle from the Algoim library. The hyperrectangle is templated on the type real and the dimension dim.

Parameters:

rectangle – A constant reference to an algoim::hyperrectangle object representing the bounding box dimensions and coordinates.

void set(const Point<dim> &min, const Point<dim> &max)

Sets the bounds of the domain.

Parameters:
  • min – Minimum coordinates of the box.

  • max – Maximum coordinates of the box.

void extend(const Point<dim> &point)

Expands the current bounding box such that it contains the given point.

Parameters:

point – Point that the expanded box must contain.

Query methods

real min(int dir) const

Gets the minimum value of the box along direction dir.

Parameters:

dir – Direction along which the minimum value is retrieved.

Returns:

Minimum value along dir.

real max(int dir) const

Gets the maximum value of the box along direction dir.

Parameters:

dir – Direction along which the maximum value is retrieved.

Returns:

Maximum value along dir.

template<int dim_aux = dim>
real min() const

Gets the minimum value of the box for 1D boxes.

Returns:

Minimum value.

template<int dim_aux = dim>
real max() const

Gets the maxnimum value of the box for 1D boxes.

Returns:

Maximum value.

const Point<dim> &min_corner() const

Gets the minimum bounds along all the directions.

Returns:

Minimum bounds along all dim directions.

const Point<dim> &max_corner() const

Gets the maximum bounds along all the directions.

Returns:

Maximum bounds along all dim directions.

BoundBox<dim> extend(real delta) const

Extends the current bounding box by a given +/- delta on each side.

Parameters:

delta – Amount by which the box is extended.

Returns:

Extended bounding box.

::algoim::HyperRectangle<real, dim> to_hyperrectangle() const

Converts the current object to an Algoim’s hyperrectangle.

This function transforms the current object into an instance of algoim::hyperrectangle with the specified real type and dimension.

Returns:

An instance of algoim::hyperrectangle<real, dim> representing the current object as a hyperrectangle.

Point<dim> get_lengths() const

Retrieves the lengths of the bounding box along each dimension.

This function returns a Point object representing the lengths of the bounding box along each dimension.

Returns:

A Point object containing the lengths of the bounding box along each dimension.

real length(int dir) const

Gets the box length along the given direction.

Parameters:

dir – Direction along which the length is computed.

Returns:

Box’s length along dir.

template<int dim_aux = dim>
real length() const

Gets the box length for 1D boxes.

Returns:

Box’s length.

real volume() const

Computes the box volume.

Returns:

Box’s volume.

Point<dim> mid_point() const

Gets the mid point of the box.

Returns:

Box’s mid point.

Point<dim> scale_to_new_domain(const BoundBox<dim> &new_domain, const Point<dim> &point) const

Scales the given point from the current domain to the new_domain.

Parameters:
  • new_domain – New domain to which the point is scaled to.

  • point – Point to be scaled.

Returns:

Scaledd point.

Point<dim> scale_to_0_1(const Point<dim> &point) const

Scales the given point from the current domain to the [0,1]^dim domain.

Parameters:

point – Point to be scaled.

Returns:

Scaled point.

Point<dim> scale_from_0_1(const Point<dim> &point_01) const

Scales the given point from the [0,1]^dim domain to the current domain.

Parameters:

point_01 – Point to be scaled.

Returns:

Scaled point.

template<int dim_aux = dim>
inline BoundBox<dim - 1> slice(const int const_dir) const

Performs a slice of the domain reduding it by one dimension.

Returns:

Dim-1 dimensional domain.

template<int dim>
class CartGridTP
#include <cart_grid_tp.hpp>

Class representing a dim-dimensional Cartesian tensor-product grid.

Template Parameters:

dim – Dimension of the grid’s domain.

Constructors

explicit CartGridTP(const std::array<std::vector<real>, dim> &breaks)

Construct a new CartGridTP object from its breaks.

Parameters:

breaks – Breaks defining the cell intervals along the dim parametric directions.

CartGridTP(const BoundBox<dim> &domain, const std::array<std::size_t, dim> &n_intrvs_dir)

Construct a new CartGridTP object from a domain and the number of intervals per direction.

Parameters:
  • domain – Domain of the grid.

  • n_intrvs_dir – Number of cell intervals along the dim parametric directions.

explicit CartGridTP(const std::array<std::size_t, dim> &n_intrvs_dir)

Construct a new CartGridTP object from a [0,1] domain and the number of intervals per direction.

Parameters:

n_intrvs_dir – Number of cell intervals along the dim parametric directions.

Query methods

const std::vector<real> &get_breaks(int dir) const

Gets the breaks along the given direction dir.

Parameters:

dir – Direction along which the breaks are extracted.

Returns:

Breaks along dir.

const BoundBox<dim> &get_domain() const

Gets the grid’s domain.

Returns:

Grid’s domain.

int get_cell_id(const PointType &point, const Tolerance &tolerance = Tolerance()) const

Get the id of the cell the given point belongs to.

Parameters:
  • point – Point to query.

  • toleranceTolerance to be used in checkings.

Returns:

Id of the cell. If the point belongs to more than one cell, it returns the lowest index of the neighbor cells.

std::optional<int> at_cells_boundary(const PointType &point, const Tolerance &tolerance = Tolerance()) const

Checks if the given point is on the boundary of two cells (up to tolerance).

Parameters:
  • point – Point to query.

  • toleranceTolerance to be used in checkings.

Returns:

If the point is not at any boundary, the returned optiona type has no value. Otherwise, it contains the index of the constant direction (from 0 to dim-1) of the boundary. For instance, in a 2D grid, if the point lays in a vertical boundary between two cells, it returns 0, but if the boundary is horizontal, returns 1.

bool on_boundary(int cell_id, int local_facet_id) const

Checks if a cell’s facet is on the grids boundary.

Parameters:
  • cell_id – Id of the cell whose facet is checked.

  • local_facet_id – Id of the local facet of the cell.

Returns:

bool Whether the face is on the grid’s boundary.

std::vector<int> get_boundary_cells(int facet_id) const

Gets the list of cells belonging to given facet of the grid.

Parameters:

facet_id – Id of the grid facet. It must be in the range [0, dim*2).

Returns:

List of cells on the facet.

int to_flat(const TensorIndexTP<dim> &tid) const

Gets the flat index of a grid cell from the tensor index.

Parameters:

tid – Tensor cell index to transform.

Returns:

Flat cell index.

TensorIndexTP<dim> to_tensor(int fid) const

Gets the tensor index of a grid cell from the flat index.

Parameters:

fid – Flat cell index to transform.

Returns:

Tensor cell index.

TensorSizeTP<dim> get_num_cells_dir() const

Gets the number of cells per direction.

Returns:

Number of cells per direction.

int get_num_cells() const

Gets the total number of cell.

Returns:

Total number of cells.

BoundBox<dim> get_cell_domain(int cell_fid) const

Gets an cell’s domain.

Parameters:

cell_fid – Flat id of the cell.

Returns:

Bounding box of the cell’s domain.

template<int dim>
struct CutCellsQuad
#include <cut_quadrature.hpp>

Public Members

std::vector<int> cells

List of cut cells.

std::vector<int> n_pts_per_cell

List of number of quadrature points per cut cell.

std::vector<Point<dim>> points

Vector of points.

std::vector<real> weights

Vector of weights.

template<int dim>
struct CutIsoBoundsQuad
#include <cut_quadrature.hpp>

Public Members

std::vector<int> cells

List of cut cells.

std::vector<int> local_facet_ids

List of (local) facet ids associated to the vector cells.

std::vector<int> n_pts_per_facet

List of number of quadrature points per cut facet.

std::vector<Point<dim>> points

Vector of points.

std::vector<real> weights

Vector of weights.

template<int dim>
struct CutUnfBoundsQuad
#include <cut_quadrature.hpp>

Public Members

std::vector<int> cells

List of cut cells.

std::vector<int> n_pts_per_cell

List of number of quadrature points per cut cell.

std::vector<Point<dim>> points

Vector of points.

std::vector<real> weights

Vector of weights.

std::vector<Point<dim>> normals

Vector of normals.

class Gauss
#include <quadrature.hpp>

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.

Public Static Functions

static real get_abscissa(const int n_points, const int pt_id)

Gets a Gauss-Legendre abscissa in the [0,1] domain.

Parameters:
  • n_points – Number of points of the 1D quadrature. It must be smaller or equal than n_max.

  • pt_id – Id of the point whose abscissa is returned. The point must be in the range [0, n_points [.

Returns:

Abscissa of the point.

static real get_weight(const int n_points, const int pt_id)

Gets a Gauss-Legendre weight in the [0,1] domain.

Parameters:
  • n_points – Number of points of the 1D quadrature. It must be smaller or equal than n_max.

  • pt_id – Id of the point whose weight is returned. The point must be in the range [0, n_points [.

Returns:

Weight of the point.

Public Static Attributes

static int n_max = 100

Maximum precomputed degree of Gauss-Legendre points.

template<int dim>
class Quadrature
#include <quadrature.hpp>

Class for storing dim-dimensional quadratures (non-tensor product).

Public Types

using Pt = Point<dim>

Point type.

using Domain = std::conditional_t<dim == 1, Interval, BoundBox<dim>>

Domain type.

Public Functions

Quadrature() = delete

Default constructor.

Warning

Not allowed to be used.

Quadrature(const std::vector<Pt> &points, const std::vector<real> &weights)

Constructor.

Note

points and weight must have the same length.

Parameters:
  • points – Vector of points referred to [0, 1].

  • weights – Vector of weights referred to [0, 1].

void scale_weights(const real ratio)

Scales the weights multiplying them by the given ratio.

Parameters:

ratio – Ratio respect to which the weights are scaled.

const std::vector<Pt> &points() const

Returns the vector of points.

Returns:

Constant reference to the vector of points.

const std::vector<real> &weights() const

Returns the vector of weights.

Returns:

Constant reference to the vector of weights.

std::vector<Pt> &points()

Returns the vector of points.

Returns:

Non-constant reference to the vector of points.

std::vector<real> &weights()

Returns the vector of weights.

Returns:

Non-constant reference to the vector of weights.

std::size_t get_num_points() const

Returns the number of points (and weights).

Returns:

Number of points.

Public Static Functions

static std::shared_ptr<Quadrature> create_Gauss_01(int n_points)

Creates a new quadrature class instance with a Gauss-Legendre quadrature in [0, 1].

Parameters:

n_points – Number of points of the quadrature rule along all the parametric directions.

Returns:

Gauss-Legendre quadrature rule wrapped in a shared pointer.

static std::shared_ptr<Quadrature> create_Gauss_01(const std::array<int, dim> n_points)

Creates a new quadrature class instance with a Gauss-Legendre quadrature in [0, 1].

Parameters:

n_points – Number of points of the quadrature rule in each parametric direction.

Returns:

Gauss-Legendre quadrature rule wrapped in a shared pointer.

static std::shared_ptr<Quadrature> create_Gauss(const std::array<int, dim> n_points, const Domain &domain)

Creates a new quadrature class instance with a Gauss-Legendre quadrature in the domain box.

Parameters:
  • n_points – Number of points of the quadrature rule in each parametric direction.

  • domain – Domain in which the quadrature is created.

Returns:

Gauss-Legendre quadrature rule wrapped in a shared pointer.

static std::shared_ptr<Quadrature> create_tanh_sinh_01(int n_points)

Creates a new quadrature class instance with a tanh-sinh quadrature in [0, 1].

Parameters:

n_points – Number of points of the quadrature rule along all the parametric directions.

Returns:

Tanh-sinh quadrature rule wrapped in a shared pointer.

static std::shared_ptr<Quadrature> create_tanh_sinh_01(const std::array<int, dim> n_points)

Creates a new quadrature class instance with a tanh-sinh quadrature in [0, 1].

Parameters:

n_points – Number of points of the quadrature rule in each parametric direction.

Returns:

Tanh-sinh quadrature rule wrapped in a shared pointer.

static std::shared_ptr<Quadrature> create_tanh_sinh(const std::array<int, dim> n_points, const Domain &domain)

Creates a new quadrature class instance with a tanh-sinh quadrature in the domain box.

Parameters:
  • n_points – Number of points of the quadrature rule in each parametric direction.

  • domain – Domain in which the quadrature is created.

Returns:

Tanh-sinh quadrature rule wrapped in a shared pointer.

template<int dim, int range>
class ReparamMesh
#include <reparam_mesh.hpp>

Class for storing an implicit domain reparameterization using Lagrange cells.

This reparameterization is non-conforming (it may contain hanging nodes), but is intended to be used for visualization purposes.

It stores the cells as a list of point and the connectivity. It also stores the connectivity for the wirebasket of the reparameterization.

Template Parameters:
  • dim – Parametric dimension of the reparameterization.

  • range – Range (or physical dimension) of the reparameterization.

Subclassed by qugar::impl::ImplReparamMesh< dim, range >

Public Functions

explicit ReparamMesh(int order)

Constructor.

Parameters:

order – Reparameterization order (number of points per direction).

virtual ~ReparamMesh() = default

Default virtual destructor.

void merge(const ReparamMesh<dim, range> &mesh)

Merges the given mesh into the current mesh.

This function takes another mesh as input and merges it with the current mesh. The points, connectivity, and wires connectivity are appended, shifting the indices by the number of points in the current mesh.

Note

Coincident points are not merged and duplicated wires are not purged.

Warning

Both meshes must have the same order.

Parameters:

mesh – The mesh to be merged with the current mesh.

template<int aux_dim = dim>
void add_full_cell(const BoundBox<dim> &domain, bool wirebasket)

Adds a full cell to the mesh corresponding to the given domain.

Parameters:
  • domain – The bounding box that defines the domain to which the full cell will be added.

  • wirebasket – Whether wirebasket for each cell should be added.

template<int aux_dim = dim>
void add_full_cells(const CartGridTP<dim> &grid, const std::vector<int> &cell_ids, bool wirebasket)

Adds full cells to the mesh corresponding to the grid.

Note

No coincidents points are merged. If needed, merge_coincident_points should be called after this method.

Parameters:
  • grid – Cartesian grid that define the cells to be added.

  • cell_ids – List of cells to add.

  • wirebasket – Whether wirebasket for each cell should be added.

bool use_Chebyshev() const

Determines whether the Chebyshev nodes are used, or equally spaced nodes.

Chebyshev nodes should be used for high-order reparameterizations if the reparameterization order is greater or equal than chebyshev_order.

Returns:

true if the Chebyshev method is used, false otherwise.

void insert_cell_point(const Point<range> &point, int cell_id, int pt_id)

This method sets the given point point, with index pt_id, to the cell designed by cell_id.

Parameters:
  • point – Point to be set.

  • cell_id – Flat index of the cell in which the point is set.

  • pt_id – Flat index of the reparameterization point referred to the cell.

void merge_coincident_points(const Tolerance &tol)

Merges coincident points in the reparameterization up to tolerance.

Note

The connecitivity of cells (and wirebasket) is updated accordingly, and duplicated wirebasket edges are purged.

Parameters:

tolTolerance to be used in the comparisons between points.

void scale_points(const std::vector<int> &point_ids, const BoundBox<range> &old_domain, const BoundBox<range> &new_domain)

Scales points from an old domain to a new domain.

This function takes two bounding boxes representing the old and new domains, and scales the points from the old to the new one.

Parameters:
  • point_ids – List of points to scale.

  • old_domain – The bounding box representing the old domain.

  • new_domain – The bounding box representing the new domain.

void scale_points(const BoundBox<range> &old_domain, const BoundBox<range> &new_domain)

Scales points from an old domain to a new domain.

This function takes two bounding boxes representing the old and new domains, and scales the points from the old to the new one.

Parameters:
  • old_domain – The bounding box representing the old domain.

  • new_domain – The bounding box representing the new domain.

const std::vector<Point<range>> &get_points() const

Retrieves the points the reparameterization points.

Returns:

A constant reference to a vector containing points.

const std::vector<std::size_t> &get_connectivity() const

Retrieves the reparameterization connectivity.

Returns:

A constant reference to the connectivity vector.

const std::vector<std::size_t> &get_wires_connectivity() const

Retrieves the wirebasket reparameterization connectivity.

Returns:

A constant reference to the wirebasket connectivity vector.

int get_order() const

Retrieves the reparameterization’s order.

Returns:

int Reparameterization order.

std::size_t get_num_cells() const

Retrieves the number of reparameterization cells.

Returns:

The number of reparameterization cells.

std::size_t get_num_points() const

Retrieves the number of points.

Returns:

The number of points.

std::size_t get_num_points_per_cell() const

Retrieves the number of points per cell, that depends on the reparameterization order.

Returns:

The number of points per cell as a std::size_t.

void write_VTK_file(const std::string &filename) const

Writes the reparameterization data to a VTK file.

This function takes a reparameterization object and writes its data to a file in VTK format, which can be used for visualization in tools like ParaView.

If the repameterization contains a wirebasket mesh, the generated VTK files is composed of three files: a file with extension .vtmb that calls two .vtu files (one for the internal reparameterization), and one for the wirebasket.

If no wirebasket mesh is present, the generated VTK file is a single .vtu file for the internal reparameterization.

The reparameterizations of both the internal and wirebasket meshes is written as high-order Lagrange VTK cells.

Parameters:

filename – The name of the file to which the data will be written, without extension.

void permute_cell_directions(std::size_t cell_id)

Permutes the directions of the given cell.

If dimension is greater or equal than 2, the first two direction are permuted. Otherwise, if the dimension is 1, the first direction is reversed.

Parameters:

cell_id – The identifier of the cell whose directions are to be permuted.

void reserve_cells(std::size_t n_new_cells)

Reserves memory for a specified number of cells.

It only reserves memory, no resize performed, for poitns and cells connectivity. Not for wires.

Parameters:

n_new_cells – The number of cells to allocate memory for.

int allocate_cells(std::size_t n_new_cells)

Allocates memory for a specified number of cells.

For the points, it just reserves memory for the new points. Thus, new points should be pushed back.

However, for the cells connectivity, it resizes the array, ready to insert new indices at the corresponding positions.

Parameters:

n_new_cells – The number of cells to allocate memory for.

Returns:

Id of the first allocated cell.

Public Static Functions

static bool use_Chebyshev(int order)

Determines whether the Chebyshev nodes are used, or equally spaced nodes.

Chebyshev nodes should be used for high-order reparameterizations if the reparameterization order is greater or equal than chebyshev_order.

Returns:

order Reparameterization order.

Returns:

true if the Chebyshev method is used, false otherwise.

Public Static Attributes

static const int chebyshev_order = 7

Order for deciding if Chebyshev nodes are used for defining the Lagrange polynomials instead of equally spaced nodes. If the reparameterization order is greater than this value, Chebyshev nodes are used.

template<int dim>
class SubCartGridTP
#include <cart_grid_tp.hpp>

Subgrid of a Cartesian grid TP. It is a subset of the cells of a given grid.

Template Parameters:

dim – Parametric dimension.

Public Functions

SubCartGridTP(const CartGridTP<dim> &grid, const TensorIndexTP<dim> &indices_start, const TensorIndexTP<dim> &indices_end)

Constructor.

Parameters:
  • grid – Parent grid.

  • indices_start – Start indices of the coordinates of the parent grid.

  • indices_end – End indices of the coordinates of the parent grid.

SubCartGridTP(const CartGridTP<dim> &grid, const TensorIndexRangeTP<dim> &indices_range)

Constructor.

Parameters:
  • grid – Parent grid.

  • indices_range – Indices range.

explicit SubCartGridTP(const CartGridTP<dim> &grid)

Constructor. Creates a subgrid containing the full grid.

Parameters:

grid – Parent grid.

TensorSizeTP<dim> get_num_cells_dir() const

Gets the number of cells (spans) per direction of the subgrid.

The ordering is such that dimension 0 is inner-most, i.e., iterates the fastest, while dimension dim-1 is outer-most and iterates the slowest.

Returns:

Number of cells per direction.

int to_flat(const TensorIndexTP<dim> &tid) const

Gets the flat index of a grid cells from the tensor index.

Parameters:

tid – Tensor cell index to transform.

Returns:

Flat cell index.

int get_num_cells() const

Gets the total number of cells of the subgrid.

Returns:

Total number of cells.

bool is_unique_cell() const

Checks if the subgrid has only one cell.

Returns:

True if it has only one cell, false otherwise.

const TensorIndexRangeTP<dim> &get_range() const

Gets a range describing the range of the subrid.

Returns:

Tensor-product indices range.

std::array<std::shared_ptr<const SubCartGridTP<dim>>, 2> split() const

Splits the current subgrid along the direction with a largest number of cells.

Returns:

Two generated subgrid wrapped in shared pointers.

BoundBox<dim> get_domain() const

Creates the bounding box of the subgrid’s domain.

Returns:

Bounding box of the subgrid.

const CartGridTP<dim> &get_grid() const

Gets the parent grid.

Returns:

Constant reference to the parent grid.

int get_single_cell() const

Gets the single cell in the subgrid.

Returns:

Single cell in the subgrid.

class TanhSinh
#include <quadrature.hpp>

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.

Public Static Functions

static real get_abscissa(const int n_points, const int pt_id)

Gets a tanh-sinh abscissa in the [0,1] domain.

Parameters:
  • n_points – Number of points of the 1D quadrature. It must be smaller or equal than n_max.

  • pt_id – Id of the point whose abscissa is returned. The point must be in the range [0, n_points [.

Returns:

Abscissa of the point.

static real get_weight(const int n_points, const int pt_id)

Gets a tanh-sinh weight in the [0,1] domain.

Parameters:
  • n_points – Number of points of the 1D quadrature. It must be smaller or equal than n_max.

  • pt_id – Id of the point whose weight is returned. The point must be in the range [0, n_points [.

Returns:

Weight of the point.

Public Static Attributes

static int n_max = 100

Maximum precomputed degree of tanh-sinh points.

template<int dim>
class TensorIndexRangeTP
#include <tensor_index_tp.hpp>

Class representing a dim-dimensional range defined by lower and upper tensor bounds.

Template Parameters:

dim – Dimension of the indices.

Public Functions

TensorIndexRangeTP(const TensorIndexTP<dim> &lower_bound, const TensorIndexTP<dim> &upper_bound)

Constructs a new TensorIndexRangeTP object from its lower and upper bounds.

Parameters:
  • lower_bound – Lower bounds.

  • upper_bound – Upper bounds.

explicit TensorIndexRangeTP(const TensorIndexTP<dim> &upper_bound)

Constructs a new TensorIndexRangeTP object from its upper bound. The lower bound is assumed to be zero.

Parameters:

upper_bound – Upper bounds.

explicit TensorIndexRangeTP(const TensorSizeTP<dim> &upper_bound)

Constructs a new TensorIndexRangeTP object from its upper bound. The lower bound is assumed to be zero.

Parameters:

upper_bound – Upper bounds.

explicit TensorIndexRangeTP(int upper_bound)

Constructs a new TensorIndexRangeTP object from its upper bound. The lower bound is assumed to be zero.

Parameters:

upper_bound – Upper bounds. All components are set to the given value.

TensorSizeTP<dim> get_sizes() const

Gets the sizes along all the directions.

Returns:

Range sizes along all the directions.

int size() const

returns the number of entries in the range.

Returns:

Number of entries in the range.

std::array<TensorIndexRangeTP<dim>, 2> split() const

Splits the current range along the direction with a largest number of indices.

Returns:

Two generated ranges wrapped in shared pointers.

Iterator cbegin() const

Creates a begin iterator.

Note

Constant version.

Returns:

Created iterator.

Iterator begin() const

Creates a begin iterator.

Returns:

Created iterator.

Iterator cend() const

Creates an end iterator.

Note

Constant version.

Returns:

Created iterator.

Iterator end() const

Creates an end iterator.

Returns:

Created iterator.

const TensorIndexTP<dim> &get_lower_bound() const

Gets the lower bound.

Returns:

Lower bound.

const TensorIndexTP<dim> &get_upper_bound() const

Gets the upper bound.

Returns:

Upper bound.

bool is_in_range(const TensorIndexTP<dim> &index) const

Checks that the given index is contained in the range.

Parameters:

index – Index to be checked.

Returns:

Whether the index is contained in the range or not.

bool is_in_range(int index) const

Checks that the given index is contained in the range.

Parameters:

index – Index to be checked.

Returns:

Whether the index is contained in the range or not.

class Iterator
#include <tensor_index_tp.hpp>

Iterator class for tensor index ranges.

Public Functions

Iterator(const TensorIndexTP<dim> &index, const TensorIndexTP<dim> &lower_bound, const TensorIndexTP<dim> &upper_bound)

Constructs a new iterator object given and index and its lower and upper bounds.

Parameters:
  • index – Current index for the iterator.

  • lower_bound – Lower bounds.

  • upper_bound – Upper bounds.

const TensorIndexTP<dim> &operator*() const

Dereference operator.

Returns:

Reference to the index pointed by the iterator.

const TensorIndexTP<dim> *operator->() const

Indirection operator.

Returns:

Pointer to the index pointed by the iterator.

int flat() const

Transforms the current tensor index into a flat one by considering the upper bound as associated size.

Returns:

Computed flat index.

Iterator &operator++()

Pre increments the iterator.

Returns:

Updated iterator.

Iterator operator++(int)

Post increments the iterator.

Returns:

Updated iterator.

bool operator==(const Iterator &rhs) const

Compares if two iterators are equal. They are considered to be equal if their indices and bounds are equal for all the components.

Parameters:

rhsIterator to compare to the current one.

Returns:

Whether both iterators are equal.

bool operator!=(const Iterator &rhs) const

Compares if two iterators are different. They are considered to be different if at least one of the components of their indices or bounds are different.

Parameters:

rhsIterator to compare to the current one.

Returns:

Whether both iterators are different.

template<int dim>
class TensorIndexTP : public Vector<int, dim>
#include <tensor_index_tp.hpp>

Class representing a dim-dimensional tensor-product indices.

Template Parameters:

dim – Dimension of the indices.

Constructors

template<int aux_dim = dim>
inline TensorIndexTP()

Default constructor. Initializes indices to zero.

template<int aux_dim = dim>
inline explicit TensorIndexTP(int ind)

Constructs a new TensorIndexTP object from an index.

Note

Not available for 0-dimensional case.

Parameters:

ind – All the indices are set to this values.

explicit TensorIndexTP(const Vector<int, dim> &indices)

Constructs a new TensorIndexTP object from its indices.

Parameters:

indices – Indices along the parametric directions.

explicit TensorIndexTP(const TensorSizeTP<dim> &indices)

Constructs a new TensorIndexTP object from its indices.

Parameters:

indices – Indices along the parametric directions.

template<int aux_dim = dim>
inline TensorIndexTP(int flat_index, const TensorIndexTP<dim> &size)

Constructs a new TensorIndexTP object from a flat index and an associated tensor size following the lexicographical ordering convention (i.e., lower indices run faster than higher ones).

Note

Not available for 0-dimensional case.

Parameters:
  • flat_index – Flat index.

  • size – Tensor size the flat_index is associated to.

template<int aux_dim = dim>
inline TensorIndexTP(int flat_index, const TensorSizeTP<dim> &size)

Constructs a new TensorIndexTP object from a flat index and an associated tensor size following the lexicographical ordering convention (i.e., lower indices run faster than higher ones).

Note

Not available for 0-dimensional case.

Parameters:
  • flat_index – Flat index.

  • size – Tensor size the flat_index is associated to.

template<typename ...T>
inline explicit TensorIndexTP(T... indices)

Constructs the tensor index from an argument list.

Parameters:

indices – Input arguments for the indices.

inline const Vector<int, dim> &as_Vector() const

Returns the tensor index casted as a Vector.

Note

Constant version.

Returns:

A constant reference to a Vector.

Vector<int, dim> &as_Vector()

Returns the tensor index casted as a Vector.

Note

Non-constant version.

Returns:

A constant reference to a Vector.

template<typename S>
int flat(const S &size) const

Gets the flat index associated to the current tensor index for a given tensor size.

Parameters:

size – Tensor size the current index is associated to.

Template Parameters:

S – Size type.

Returns:

Computed flat index.

bool operator==(const TensorIndexTP<dim> &rhs) const

Checks if two tensor indices are equal by comparing all its components. They are equal if all the components are equal.

Returns:

Whether they are equal or not.

bool operator!=(const TensorIndexTP<dim> &rhs) const

Checks if two tensor indices are different by comparing all its components. They are different if at least one component is different.

Returns:

Whether they are different or not.

bool operator<(const TensorIndexTP<dim> &rhs) const

Checks if the current tensor index is smaller than rhs. To determine if one is smaller than the other, all the componets are compared starting from the last one.

Returns:

Whether the current tensor is smaller than rhs.

template<int aux_dim = dim>
TensorIndexTP<dim - 1> remove_component(int comp) const

Creates a new index by removing one of its component.

Parameters:

comp – Index of the component to be removed.

Returns:

New generated index.

std::size_t hash() const

Computes a hash value for the tensor index.

This function generates a hash value that uniquely identifies the tensor index. The hash value can be used for efficient lookups and comparisons.

Returns:

A std::size_t value representing the hash of the tensor index.

template<int dim>
class TensorSizeTP : public Vector<int, dim>
#include <tensor_index_tp.hpp>

Class representing a dim-dimensional tensor-product sizes container.

Template Parameters:

dim – Dimension of the container.

Constructors

TensorSizeTP()

Default constructor. Initializes indices to zero.

explicit TensorSizeTP(const int size)

Construct a new TensorSizeTP object from size along directions.

Parameters:

size – Sizs along all the parametric directions.

explicit TensorSizeTP(const Vector<int, dim> &sizes)

Construct a new TensorSizeTP object from sizes along directions.

Parameters:

sizes – Sizes along the parametric directions.

explicit TensorSizeTP(const TensorIndexTP<dim> &sizes)

Construct a new TensorSizeTP object from sizes along directions.

Parameters:

sizes – Sizes along the parametric directions.

const Vector<int, dim> &as_Vector() const

Returns the tensor size casted as a Vector.

Note

Constant version.

Returns:

A constant reference to a Vector.

Vector<int, dim> &as_Vector()

Returns the tensor size casted as a Vector.

Note

Non-constant version.

Returns:

A constant reference to a Vector.

int size() const

Gets the total number of entries (product of sizes along all directions).

Returns:

Total number of entries.

bool operator==(const TensorSizeTP<dim> &rhs) const

Checks if two tensor sizes are equal by comparing all its components. They are equal if all the components are equal.

Returns:

Whether they are equal or not.

TensorSizeTP<dim> operator+(const TensorSizeTP<dim> &rhs) const

Adds two TensorSizeTP objects.

This operator performs element-wise addition of two TensorSizeTP objects.

Parameters:

rhs – The right-hand side TensorSizeTP object to be added.

Returns:

A new TensorSizeTP object that is the result of the addition.

class Tolerance
#include <tolerance.hpp>

Class for tolerance related computations.

Public Functions

explicit Tolerance(const real value)

Constructor.

Parameters:

valueTolerance value (it must be greater than zero).

explicit Tolerance(const real tol_0, const real tol_1)

Constructor. Creates a new class instance using the largest of the two given tolerances.

Parameters:
  • tol_0 – First tolerance.

  • tol_1 – Second tolerance.

explicit Tolerance(const Tolerance &tol_0, const Tolerance &tol_1)

Constructor. Creates a new class instance using the largest of the two given tolerances.

Parameters:
  • tol_0 – First tolerance.

  • tol_1 – Second tolerance.

explicit Tolerance()

Default constructor. It initalizes the class instance with a near epsilon tolerance.

void update(const real tol)

Resets the tolerance value as the maximum between the current tolerance_ and tol.

Parameters:

tolTolerance value to be reset.

void update(const Tolerance &tol)

Resets the tolerance value as the maximum between the current tolerance_ and tol.

Parameters:

tolTolerance value to be reset.

real value() const

Returns the real tolerance value.

Returns:

Tolerance real value.

bool is_zero(const real val) const

Check if val is zero up to tolerance.

It checks \(|val| \leq tolerance\)

Parameters:

val – Value to be checked.

Returns:

Whether val is equal to zero up to tolerance.

bool is_negative(const real val) const

Check if val is negative.

It checks \(val \leq -tolerance\)

Parameters:

val – Value to be checked.

Returns:

Whether val is negative up to tolerance.

bool is_positive(const real val) const

Check if val is positive.

It checks \(tolerance < val\)

Parameters:

val – Value to be checked.

Returns:

Whether val is positive up to tolerance.

bool equal(const real lhs, const real rhs) const

Compares if two values are equal up to tolerance.

It checks \(|lhs -rhs| \leq tolerance\)

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether both arguments are equal up to tolerance.

bool equal_rel(const real lhs, const real rhs) const

Compares if two values are equal up to tolerance relative to the larger of the two arguments.

It checks \(|lhs -rhs| \leq tolerance \max (|lhs|, |rhs|) \)

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether both arguments are relatively equal up to tolerance.

bool equal_rel(const real lhs, const real rhs, const Tolerance &rel_tolerance) const

Compares if two values are equal up to the different absolute and relative tolerances.

It checks \(|lhs -rhs| \leq (tolerance + rel_tolerance \max (|lhs|, |rhs|)) \) where \(tolerance\) is the current tolerance and \(rel_tolerance\) is the one provided as input.

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

  • rel_tolerance – Relative tolerance.

Returns:

Whether both arguments are relatively equal up to tolerance.

bool greater_than(const real lhs, const real rhs) const

Compares if a value is greater than other up to tolerance.

It checks \((lhs - rhs) > tolerance\)

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether lhs is greater than rhs up to tolerance.

bool greater_equal_than(const real lhs, const real rhs) const

Compares if a value is greater or equal than other up to tolerance.

It checks \((lhs - rhs) > tolerance\) or \(|rhs - lhs| < tolerance\).

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether lhs is greater or equal than rhs up to tolerance.

bool greater_than_rel(const real lhs, const real rhs) const

Compares if a value is greater than other up to tolerance relative to the larger of the two argumetns.

It checks \((lhs - rhs) > tolerance \max (|lhs|, |rhs|)\)

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether lhs is relatively greater than rhs up to tolerance.

bool smaller_than(const real lhs, const real rhs) const

Compares if a value is smaller than other up to tolerance.

It checks \((rhs - lhs) > tolerance\)

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether lhs is smaller than rhs up to tolerance.

bool smaller_equal_than(const real lhs, const real rhs) const

Compares if a value is smaller or equal than other up to tolerance.

It checks \((rhs - lhs) > tolerance\) or \(|rhs - lhs| < tolerance\).

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether lhs is smaller or equal than rhs up to tolerance.

bool smaller_than_rel(const real lhs, const real rhs) const

Compares if a value is smaller than other up to tolerance relative to the larger of the two argumetns.

It checks \((rhs - lhs) > tolerance \max (|lhs|, |rhs|)\)

Parameters:
  • lhs – Left-hand-side value in the comparison.

  • rhs – Right-hand-side value in the comparison.

Returns:

Whether lhs is relatively smaller than rhs up to tolerance.

void unique(std::vector<real> &values) const

Makes the given values unique up to tolerance. values will be also sorted.

Parameters:

values – Vector of value of make unique.

template<int dim, typename T>
inline bool coincident(const Point<dim, T> &pt_0, const Point<dim, T> &pt_1) const

Checks if two points are coincident.

This function compares two points of the same dimension and type to determine if they are coincident. Two points are considered coincident if all their corresponding coordinates are equal up to tolerance.

Template Parameters:
  • dim – The dimension of the points.

  • T – The type of the coordinates of the points.

Parameters:
  • pt_0 – The first point to compare.

  • pt_1 – The second point to compare.

Returns:

true if the points are coincident, false otherwise.

template<int dim>
class UnfittedDomain
#include <unfitted_domain.hpp>

Subclassed by qugar::impl::UnfittedImplDomain< dim >

namespace impl

Typedefs

template<int dim>
using ScalarFunc = DomainFunc<dim, 1>

Alias for scalar functions.

template<int dim>
using ImplicitFunc = ScalarFunc<dim>

Alias for implicit functions.

Enums

enum FuncSign

Values:

enumerator negative
enumerator positive
enumerator 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.

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:
  • dim – The dimension of the BezierTP object.

  • range – The range of the BezierTP object.

Parameters:

func – The DomainFunc object to be checked.

Returns:

true if the object is of type BezierTP, false otherwise.

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 Parameters:

T – Type of the input coordinate.

Parameters:
  • point – Evaluation point.

  • order – Order of the polynomial.

  • values – Computed basis values.

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 Parameters:
  • T – Type of the input coordinate.

  • V – Type of the output vector of basis values.

Parameters:
  • point – Evaluation point.

  • order – Order of the polynomial.

  • der – Order of the derivative to be computed. If 0, the value itself is computed.

  • values – Computed basis values.

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.

Parameters:
  • lhs – First Bezier to multiply.

  • rhs – Second Bezier to multiply.

Returns:

Product of 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 Parameters:
  • dim – Dimension of the resultant Bezier.

  • range – Range of the resultant Bezier.

  • dim2 – Common dimension between Beziers.

Parameters:
  • lhs – Bezier to be composed.

  • rhs – Composing Bezier.

Returns:

Resultant composition.

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.

Note

The generated reparameterization has a wirebasket associated, but coincident were not merged.

Template Parameters:
  • dim – Parametric dimension of the function.

  • S – Flag 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:
  • bzr – Bezier implicit function to reparameterize.

  • domain – Domain to which the implicit function refers to (even if Beziers are defined in the unit domain).

  • order – Order of the reparameterization (number of points per direction in each reparameterization cell).

Returns:

Generated reparameterization.

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.

Note

The generated reparameterization has a wirebasket associated, but coincident were not merged.

Template Parameters:
  • dim – Parametric dimension of the function.

  • S – Flag 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:
  • bzr – Bezier implicit function to reparameterize.

  • domain – Domain to which the implicit function refers to (even if Beziers are defined in the unit domain).

  • reparam – Reparameterization container to which new generated cells are appended to.

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.

Note

The generated reparameterization has a wirebasket associated, but coincident were not merged.

Template Parameters:
  • dim – Parametric dimension of the function.

  • S – Flag 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:
  • func – Function to reparameterize.

  • domain – Domain to reparameterize.

  • order – Order of the reparameterization (number of points per direction in each reparameterization cell).

Returns:

Generated reparameterization.

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.

Note

The generated reparameterization has a wirebasket associated, but coincident were not merged.

Template Parameters:
  • dim – Parametric dimension of the function.

  • S – Flag 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:
  • func – Function to reparameterize.

  • domain – Domain to reparameterize.

  • reparam – Reparameterization container to which new generated cells are appended to.

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.

Note

The generated reparameterization has a wirebasket associated, but coincident were not merged.

Template Parameters:

dim – Parametric dimension of the function.

Parameters:
  • func – Function to reparameterize.

  • domain – Domain to reparameterize.

  • facet_id – Id of the face to reparameterize. It must be a value in the range [0, 2*dim[.

  • order – Order of the reparameterization (number of points per direction in each reparameterization cell).

Returns:

Generated reparameterization.

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.

Note

The generated reparameterization has a wirebasket associated, but coincident were not merged.

Template Parameters:

dim – Parametric dimension of the function.

Parameters:
  • func – Function to reparameterize.

  • domain – Domain to reparameterize.

  • facet_id – Id of the face to reparameterize. It must be a value in the range [0, 2*dim[.

  • reparam – Reparameterization container to which new generated cells are appended to.

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 Parameters:

dim – Dimension of the function and point.

Parameters:
  • phi – Implicit function to be tested.

  • point – Point 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.

template<int dim>
int get_facet_constant_dir(int local_facet_id)

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_id – Id 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[.

template<int dim>
int get_facet_side(int local_facet_id)

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_id – Id 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.

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.

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:

dim – The dimension of the space.

Parameters:
  • const_dir – The constant direction for which the local facet ID is to be determined.

  • side – The side (0 or 1) in the specified direction.

Returns:

The local facet ID corresponding to the given 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.

The edges are assumed to follow a lexicographical numbering.

Parameters:

local_edge_id – Id of the edge referred to a cell.

Returns:

Constant directions in the range [0,dim[.

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.

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_id – Id of the edge referred to a cell.

Returns:

Sides of the edge.

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.

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:
  • point – The point at which to evaluate the Lagrange basis polynomials.

  • order – The order of the Lagrange basis polynomials (degree + 1).

  • chebyshev – A boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).

  • values – A 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.

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.

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:
  • point – The point at which to evaluate the derivative.

  • order – The order of the Lagrange basis polynomial (degree + 1).

  • chebyshev – A boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).

  • values – A 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.

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.

Parameters:
  • point – The point at which to evaluate the Lagrange basis functions.

  • order – The order (degree + 1) of the tensor product along each direction.

  • chebyshev – A boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).

  • basis – A 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.

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.

Parameters:
  • point – The point at which to evaluate the derivative of the Lagrange basis functions.

  • order – The order (degree + 1) of the tensor product along each direction.

  • chebyshev – A boolean flag indicating whether to use Chebyshev nodes (true) of 2nd kind or standard equidistant nodes (false).

  • basis_ders – A 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.

template<int dim>
class AffineTransf
#include <affine_transf.hpp>

Class for representing affine transformations.

It transforms points as y = A * x + b, where x is the original point, y is the new point, b is the translation vector, and A is the transformation matrix that may be defined by rotations and scaling.

The member coefs_ stores the involved coefficients. Thus, first dim*dim values of coefs_ contain the values of the matrix A (row-wise), and the last dim values, the ones of the vector b.

Template Parameters:

dim – Parametric direction.

Public Functions

explicit AffineTransf()

Default constructor. Creates the identity transformation.

explicit AffineTransf(const Point<dim> &origin)

Constructs a new affine tranformation that simply translates the origin.

Parameters:

origin – New origin of the reference system.

AffineTransf(const Point<dim> &origin, real scale)

Constructs a new affine transformation that first applies an isotropic scaling and then translates the origin.

Parameters:
  • origin – New origin of the reference system.

  • scale – Scaling factor to be applied along all the directions.

template<int dim_aux = dim>
AffineTransf(const Point<dim> &origin, const Point<dim> &axis_x)

Constructs a new 2D affine transformation defined by an origin and the axis_x.

It first rotates the axes and then performs the translation.

The y-axis is computed by rotating axis_x 90 degrees counter-clockwise.

Parameters:
  • origin – New origin of the reference system.

  • axis_x – New x-axis of the system (it will be normalized).

template<int dim_aux = dim>
AffineTransf(const Point<dim> &origin, const Point<dim> &axis_x, const Point<dim> &axis_y)

Constructs a new 3D affine transformation defined by an origin and a couple of directions that define the xy-plane.

It first rotates the axes and then performs the translation.

The given axis_x and axis_y define the xy-plane, but they may not be orthonormal vectors. Thus, the way in which the new system directions are computed is as follows: First, the axis_x is normalized. Then, the z-axis is computed as the normalized cross-product of the normalized vector axis_x and axis_y. Then, the y-axis is recomputed as the cross-product between normalized vectors along the z- and the x-axes.

Parameters:
  • origin – New origin of the reference system.

  • axis_x – New x-axis of the system (it will be rescaled).

  • axis_y – Axis defining the xy-plane together with axis_x. It must be not parallel to axis_x.

template<int dim_aux = dim>
AffineTransf(const Point<dim> &origin, const Point<dim> &axis_x, real scale_x, real scale_y)

Constructs a new 2D affine transformation defined by an origin, the axis_x, and an orthotropic scaling along the x- and y-axes.

It first applies the scaling, then rotates the axes, and finally performs the translation.

The y-axis is computed by rotating axis_x 90 degrees counter-clockwise.

Parameters:
  • origin – New origin of the reference system.

  • axis_x – New x-axis of the system (it will be rescaled).

  • scale_x – Scaling along the new x-axis.

  • scale_y – Scaling along the new y-axis.

template<int dim_aux = dim>
AffineTransf(const Point<dim> &origin, const Point<dim> &axis_x, const Point<dim> &axis_y, real scale_x, real scale_y, real scale_z)

Constructs a new 3D affine transformation defined by an origin, a couple of directions that define the xy-plane, and an orthotropic scaling along the x-, y, and z-axes.

It first applies the scaling, then rotates the axes, and finally performs the translation.

The given axis_x and axis_y define the xy-plane, but they may not be orthonormal vectors. Thus, the way in which the new system directions are computed is as follows: First, the axis_x is normalized. Then, the z-axis is computed as the normalized cross-product of the normalized vector axis_x and axis_y. Then, the y-axis is recomputed as the cross-product between normalized vectors along the z- and the x-axes.

Parameters:
  • origin – New origin of the reference system.

  • axis_x – New x-axis of the system (it will be rescaled).

  • axis_y – Axis defining the xy-plane together with axis_x. It must be not parallel to axis_x.

  • scale_x – Scaling along the new x-axis.

  • scale_y – Scaling along the new y-axis.

  • scale_z – Scaling along the new z-axis.

explicit AffineTransf(const Vector<real, n_coefs> &coefs)

Constructs a new object through its coefficients. a new Affine Transf object.

See the class documentation for information about the coefficients.

Parameters:

coefs – Transformation coefficients.

AffineTransf<dim> operator*(const AffineTransf<dim> &rhs) const

Concatenates two affine transformation to generate a new one.

The new generated transformation is equivalent to, first, applying the given rhs transformation to a point (gradient, hessian), and then, the current one.

Parameters:

rhs – First transformation to apply.

Returns:

Concatenated transformations.

AffineTransf<dim> inverse() const

Inverts the current transformation.

Returns:

Inverted transformation.

template<typename T>
Vector<T, dim> transform_point(const Vector<T, dim> &point) const

Transforms a point from the original to the new reference system.

Template Parameters:

T – Type of the point.

Parameters:

point – Point to be transformed.

Returns:

Transformed point.

template<typename T>
Vector<T, dim> transform_vector(const Vector<T, dim> &vector) const

Transforms a vector from the original to the new reference system (without translation).

Template Parameters:

T – Type of the vector.

Parameters:

vector – Vector to be transformed.

Returns:

Transformed vector.

template<typename T>
Tensor<T> transform_tensor(const Tensor<T> &tensor) const

Transforms a (second-order symmetric) tensor from the original to the new reference system (without translation).

Template Parameters:

T – Type of the tensor.

Parameters:

tensor – Second-order symmetric tensor to be transformed.

Returns:

Transformed tnesor.

template<int dim, int range>
class BezierTP : public qugar::impl::PolynomialTP<dim, range>
#include <bezier_tp.hpp>

dim-dimensional tensor-product Bezier polynomial function.

Template Parameters:
  • dim – Dimension of the parametric domain.

  • range – Dimension of the image.

Public Types

using Parent = PolynomialTP<dim, range>

Parent type.

using CoefsType = typename Parent::CoefsType

Coefs type.

template<typename T>
using Value = typename Parent::template Value<T>

Value type.

template<typename T>
using Gradient = typename Parent::template Gradient<T>

Gradient type.

template<typename T>
using Hessian = typename Parent::template Hessian<T>

Hessian type.

template<int N>
using Interval = ::algoim::Interval<N>

Algoim’s interval alias.

Public Functions

explicit BezierTP(const TensorSizeTP<dim> &order)

Constructor.

The coefficients vector is allocated by the constructor.

Parameters:

order – Order of the polynomial along each parametric direction.

BezierTP(const TensorSizeTP<dim> &order, const CoefsType &value)

Constructs a constant value Bezier tensor product.

Parameters:
  • order – Order of the polynomial along each parametric direction.

  • value – The value to set for all the coefficients.

BezierTP(const TensorSizeTP<dim> &order, const std::vector<CoefsType> &coefs)

Constructor.

Parameters:
  • coefs – Coefficients of the polynomial. The ordering is such that dimension 0 is inner-most, i.e., iterates the fastest, while dimension dim-1 is outer-most and iterates the slowest.

  • order – Order of the polynomial along each parametric direction.

BezierTP(const BezierTP<dim, range> &bezier)

Copy constructor.

Parameters:

bezier – Bezier to copy.

explicit BezierTP(const MonomialsTP<dim, range> &monomials)

Constructor.

Constructs from monomials.

Parameters:

monomials – Monomials from which the Bezier is created.

virtual Value<real> operator()(const Point<dim> &point) const final

Evaluator operator.

Parameters:

point – Point at which the function is evaluated.

Returns:

Function value at point.

virtual Value<Interval<dim>> operator()(const Point<dim, Interval<dim>> &point) const final

Evaluator operator.

Parameters:

point – Point at which the function is evaluated.

Returns:

Function value at point.

virtual Gradient<real> grad(const Point<dim> &point) const final

Gradient evaluator operator.

Parameters:

point – Point at which the function’s gradient is evaluated.

Returns:

Function gradient at point.

virtual Gradient<Interval<dim>> grad(const Point<dim, Interval<dim>> &point) const final

Gradient evaluator operator.

Parameters:

point – Point at which the function’s gradient is evaluated.

Returns:

Function gradient at point.

virtual Hessian<real> hessian(const Point<dim> &point) const final

Hessian evaluator operator.

Parameters:

point – Point at which the function’s hessian is evaluated.

Returns:

Function hessian at point.

const ::algoim::xarray<CoefsType, dim> &get_xarray() const

Gets a constant reference to the stored xarray view of the polynomial coefficients.

Returns:

Constant view of the polynomial coefficients.

std::shared_ptr<BezierTP<dim, range>> raise_order(const TensorSizeTP<dim> &new_order) const

Raises the order of the Bezier tensor product.

This function creates a new Bezier tensor product with an increased order as specified by the new_order parameter. The new Bezier will represent the same polynomial function as the original, but with a higher order.

Parameters:

new_order – The new order for each dimension of the tensor product. Along each dimension, the new order must be greater or equal than the current order.

Returns:

A shared pointer to the new Bezier tensor product with the raised order.

std::shared_ptr<BezierTP<dim, range>> negate() const

Creates a new Bezier tensor product (TP) object that is the negation of the current object.

This function returns a shared pointer to a new BezierTP object where each control point has the same value but with negative sign.

Returns:

A shared pointer to the negated BezierTP object.

template<int range_aux = range>
void rescale_domain(const BoundBox<dim> &new_domain)

Recomputes (in-place) the Bezier coefficients for a new domain (that may not be [0, 1]).

Parameters:

new_domain – New domain for which the coefficientes are computed.

template<int range_aux = range>
FuncSign sign() const

Returns the sign of the function represented by the Bézier tensor-product using the properties of the control points convex hull.

Returns:

FuncSign The sign of the function.

std::shared_ptr<BezierTP<dim, range>> operator*(const BezierTP<dim, range> &rhs) const

Product of two Beziers.

Parameters:

rhs – Second Bezier to multiply.

Returns:

Product of Beziers.

std::shared_ptr<BezierTP<dim, range>> operator+(const BezierTP<dim, range> &rhs) const

Addition of two Beziers.

Parameters:

rhs – Second Bezier to sum.

Returns:

Sum of Beziers. The resulting Bezier has the maximum order of the two Beziers along each dimension.

std::shared_ptr<BezierTP<dim, range>> operator-(const BezierTP<dim, range> &rhs) const

Subtraction of two Beziers.

Parameters:

rhs – Bezier to subtract.

Returns:

Subtraction of Beziers. The resulting Bezier has the maximum order of the two Beziers along each dimension.

template<int sub_dim>
std::shared_ptr<BezierTP<sub_dim, range>> compose(const BezierTP<sub_dim, dim> &rhs) const

Composes the current Bezier with the given Bezier rhs.

Template Parameters:

sub_dim – Subdimension of the right-hand-size Bezier.

Parameters:

rhs – Right-hand-side Bezier for the composition.

Returns:

Generated Bezier composition.

Public Static Functions

template<typename T>
static Value<T> casteljau(const Point<dim, T> &point, typename std::vector<CoefsType>::const_iterator &coefs, const Vector<int, dim> &order)

Evaluates a Bezier polynomial using the Casteljau’s algorithm.

For the algorithm details check: https://en.wikipedia.org/wiki/De_Casteljau%27s_algorithm

The evaluation along the dimensions is performed by recursively calling this function.

Template Parameters:

T – Type of the input and output variables.

Parameters:
  • point – Evaluation point.

  • coefs – Iterator to the coefficients of the polynomial. Their ordering is the same as the coefficients members of the class.

  • order – Order of the polynomial along the dim dimensions.

Returns:

Value of the polynomial at point.

template<typename T>
static Vector<Value<T>, dim + 1> casteljau_der(const Point<dim, T> &point, typename std::vector<CoefsType>::const_iterator &coefs, const Vector<int, dim> &order)

Evaluates the gradient Bezier polynomial using the Casteljau’s algorithm.

For the algorithm details check: https://en.wikipedia.org/wiki/De_Casteljau%27s_algorithm

The evaluation along the dimensions is performed by recursively calling this function.

Template Parameters:

T – Type of the input and output variables.

Parameters:
  • point – Evaluation point.

  • coefs – Iterator to the coefficients of the polynomial. Their ordering is the same as the coefficients members of the class.

  • order – Order of the polynomial along the dim dimensions.

Returns:

The first component of the return vector corresponds to the value of the polynomial itself, while the gradient is stored in the following dim components.

template<int dim, int range>
class DomainFunc
#include <domain_function.hpp>

Domain functions.

Template Parameters:
  • dim – Parametric dimension.

  • range – Image dimension.

  • T – Type.

Subclassed by qugar::impl::PolynomialTP< dim, 1 >, qugar::impl::PolynomialTP< dim, range >, qugar::impl::funcs::AddFunctions< dim >, qugar::impl::funcs::Annulus, qugar::impl::funcs::Constant< dim >, qugar::impl::funcs::Cylinder, qugar::impl::funcs::Ellipsoid< dim >, qugar::impl::funcs::FuncWithAffineTransf< dim >, qugar::impl::funcs::Negative< dim >, qugar::impl::funcs::Plane< dim >, qugar::impl::funcs::Sphere< dim >, qugar::impl::funcs::SubtractFunctions< dim >, qugar::impl::funcs::Torus, qugar::impl::tpms::TPMSBase< dim >

Public Types

template<int N>
using Interval = ::algoim::Interval<N>

Algoim’s interval alias.

template<typename T>
using Value = std::conditional_t<range == 1, T, Vector<T, range>>

Value type.

Template Parameters:

T – Type of the input coordinates.

template<typename T>
using Gradient = Vector<Value<T>, dim>

Gradient type.

Template Parameters:

T – Type of the input coordinates.

template<typename T>
using Hessian = Vector<Value<T>, num_hessian>

Hessian (symmetric type).

Template Parameters:

T – Type of the input coordinates.

Public Functions

DomainFunc() = default

Default constructor.

DomainFunc(const DomainFunc&) = default

Default copy constructor.

DomainFunc(DomainFunc&&) = default

Default move constructor.

DomainFunc &operator=(const DomainFunc&) = default

Default copy assignment operator.

DomainFunc &operator=(DomainFunc&&) = default

Default move assignment operator.

virtual ~DomainFunc() = default

Default virtual destructor.

virtual Value<real> operator()(const Point<dim> &point) const = 0

Evaluator operator.

Note

This is a purely virtual method and must be implemented in derived classes.

Parameters:

point – Point at which the function is evaluated.

Returns:

Function value at point.

virtual Value<Interval<dim>> operator()(const Point<dim, Interval<dim>> &point) const = 0

Evaluator operator.

Note

This is a purely virtual method and must be implemented in derived classes.

Parameters:

point – Point at which the function is evaluated.

Returns:

Function value at point.

virtual Gradient<real> grad(const Point<dim> &point) const = 0

Gradient evaluator operator.

Note

This is a purely virtual method and must be implemented in derived classes.

Parameters:

point – Point at which the function’s gradient is evaluated.

Returns:

Function gradient at point.

virtual Gradient<Interval<dim>> grad(const Point<dim, Interval<dim>> &point) const = 0

Gradient evaluator operator.

Note

This is a purely virtual method and must be implemented in derived classes.

Parameters:

point – Point at which the function’s gradient is evaluated.

Returns:

Function gradient at point.

virtual Hessian<real> hessian(const Point<dim> &point) const = 0

Hessian evaluator operator.

Note

This is a purely virtual method and must be implemented in derived classes.

Parameters:

point – Point at which the function’s hessian is evaluated.

Returns:

Function hessian at point.

Public Static Attributes

static const int num_hessian = dim * (dim + 1) / 2

Number of Hessian (symmetric) components.

template<int dim, int range>
class ImplReparamMesh : public qugar::ReparamMesh<dim, range>
#include <impl_reparam_mesh.hpp>

Class for storing an implicit domain reparameterization using Lagrange cells.

This reparameterization is non-conforming (it may contain hanging nodes), but is intended to be used for visualization purposes.

It stores the cells as a list of point and the connectivity. It also stores the connectivity for the wirebasket of the reparameterization.

Template Parameters:
  • dim – Parametric dimension of the reparameterization.

  • range – Range (or physical dimension) of the reparameterization.

Public Functions

explicit ImplReparamMesh(int order)

Constructor.

Parameters:

order – Reparameterization order (number of points per direction).

void generate_wirebasket(const std::vector<std::reference_wrapper<const ImplicitFunc<range>>> &impl_funcs, const BoundBox<range> &domain, const Tolerance &tol)

Generates the wirebasket for all the cells of the current reparameterization.

Parameters:
  • impl_funcs – Implicit functions defining the domain.

  • domain – Domain (bounding box) to which the reparameterization corresponds.

  • tolTolerance to be used if points belong to the wirebasket.

void generate_wirebasket(const std::vector<std::reference_wrapper<const ImplicitFunc<range>>> &impl_funcs, const std::vector<int> &cell_ids, const BoundBox<range> &domain, const Tolerance &tol)

Generates the wirebasket for the given cell_ids of the current reparameterization.

Parameters:
  • impl_funcs – List of implicit functions defining the domain.

  • cell_ids – List of cells ids to be considered for extracting the wirebasket.

  • domain – Domain (bounding box) to which the reparameterization corresponds.

  • tolTolerance to be used if points belong to the wirebasket.

template<int dim_aux = dim>
void orient_cells_positively()

Reorients the cells in the mesh to ensure they are oriented positively.

This function modifies the orientation of the cells (by modifying its connectivity) in the mesh so that all cells have a positive orientation.

Positive orientation is considered the one such that determinant of the Jacobian (evaluated at the cell’s mid-point) is positive.

This is typically required for consistency in computational geometry and finite element methods, where the orientation of cells can affect the results of calculations.

Note

The wirebasket connectivity is not modified at all.

template<int dim_aux = dim>
void orient_cells_positively(const std::vector<int> &cell_ids)

Reorients the cells in the mesh to ensure they are oriented positively.

This function modifies the orientation of the cells (by modifying its connectivity) in the mesh so that all cells have a positive orientation.

Positive orientation is considered the one such that determinant of the Jacobian (evaluated at the cell’s mid-point) is positive.

This is typically required for consistency in computational geometry and finite element methods, where the orientation of cells can affect the results of calculations.

Note

The wirebasket connectivity is not modified at all.

Parameters:

cell_ids – List of cells ids to be considered for reorientation.

template<int dim_aux = dim>
void orient_levelset_cells_positively(const std::vector<std::reference_wrapper<const ImplicitFunc<range>>> &impl_funcs)

Reorients the levelset cells in the mesh to ensure their normal is an outer normal.

This function modifies the orientation of the levelset cells (by modifying its connectivity) in the mesh so that all cells have a positive orientation.

Positive orientation is considered the one such that the cell normal (evaluated at the cell’s mid-point) goes along the normal provided by the given implicit functions at that point. I.e., the dot product between both normals is positive.

Among all the provided functions, the one whose value at the mid-point is minimum is chosen for evaluating its normals.

This is typically required for consistency in computational geometry and finite element methods, where the orientation of cells can affect the results of calculations.

Note

The wirebasket connectivity is not modified at all.

Parameters:

impl_funcs – List of functions that define the implicit domain and according to which the orientation is computed.

template<int dim_aux = dim>
void orient_levelset_cells_positively(const std::vector<std::reference_wrapper<const ImplicitFunc<range>>> &impl_funcs, const std::vector<int> &cell_ids)

Reorients the levelset cells in the mesh to ensure their normal is an outer normal.

This function modifies the orientation of the levelset cells (by modifying its connectivity) in the mesh so that all cells have a positive orientation.

Positive orientation is considered the one such that the cell normal (evaluated at the cell’s mid-point) goes along the normal provided by the given implicit functions at that point. I.e., the dot product between both normals is positive.

Among all the provided functions, the one whose value at the mid-point is minimum is chosen for evaluating its normals.

This is typically required for consistency in computational geometry and finite element methods, where the orientation of cells can affect the results of calculations.

Note

The wirebasket connectivity is not modified at all.

Parameters:
  • impl_funcs – List of functions that define the implicit domain and according to which the orientation is computed.

  • cell_ids – List of cells ids to be considered for reorientation.

template<int dim_aux = dim>
void orient_facet_cells_positively(int local_facet_id)

Orients the cells reparameterization a domain facet positively.

This function ensures that the cells reparameterizing a facet are oriented positively according to the specified local facet ID of a hyperrectangular domain. I.e., the cells normals are oriented towards the exterior of the domain.

Parameters:

local_facet_id – The local ID of the facet to be oriented positively.

template<int dim_aux = dim>
void orient_facet_cells_positively(int local_facet_id, const std::vector<int> &cell_ids)

Orients the cells reparameterization a domain facet positively.

This function ensures that the cells reparameterizing a facet are oriented positively according to the specified local facet ID of a hyperrectangular domain. I.e., the cells normals are oriented towards the exterior of the domain.

Parameters:
  • local_facet_id – The local ID of the facet to be oriented positively.

  • cell_ids – List of cells ids to be considered for reorientation.

template<int dim, int range = 1>
class MonomialsTP : public qugar::impl::PolynomialTP<dim, 1>
#include <monomials_tp.hpp>

dim-dimensional tensor-product monomials function.

Template Parameters:
  • dim – Dimension of the parametric domain.

  • range – Dimension of the image.

Public Types

using Parent = PolynomialTP<dim, range>

Parent type.

using CoefsType = typename Parent::CoefsType

Coefs type.

template<typename T>
using Value = typename Parent::template Value<T>

Value type.

template<typename T>
using Gradient = typename Parent::template Gradient<T>

Gradient type.

template<typename T>
using Hessian = typename Parent::template Hessian<T>

Hessian type.

template<int N>
using Interval = ::algoim::Interval<N>

Algoim’s interval alias.

Public Functions

explicit MonomialsTP(const TensorSizeTP<dim> &order)

Constructor.

The coefficients vector is allocated by the constructor.

Parameters:

order – Order of the polynomial along each parametric direction.

MonomialsTP(const TensorSizeTP<dim> &order, const std::vector<CoefsType> &coefs)

Constructor.

Parameters:
  • coefs – Coefficients of the polynomial. The ordering is such that dimension 0 is inner-most, i.e., iterates the fastest, while dimension dim-1 is outer-most and iterates the slowest.

  • order – Order of the polynomial along each parametric direction.

MonomialsTP(const MonomialsTP<dim, range> &monomials)

Copy constructor.

Parameters:

monomials – Monomials to copy.

std::shared_ptr<MonomialsTP<dim, range>> create_derivative(const int dir) const

Creates a derivative of the current MonomialsTP object in the specified direction.

This function generates a new MonomialsTP object that represents the derivative of the current object in the direction specified by the parameter dir.

Note

If the monomials has already order 1 (degree 0) in the specified direction, just a clone of the current object is returned.

Parameters:

dir – The direction in which to take the derivative. This should be an integer representing the axis or dimension along which the derivative is computed.

Returns:

A std::shared_ptr to a new MonomialsTP object that represents the derivative of the current object in the specified direction.

virtual Value<real> operator()(const Point<dim> &point) const final

Evaluator operator.

Parameters:

point – Point at which the function is evaluated.

Returns:

Function value at point.

virtual Value<Interval<dim>> operator()(const Point<dim, Interval<dim>> &point) const final

Evaluator operator.

Parameters:

point – Point at which the function is evaluated.

Returns:

Function value at point.

virtual Gradient<real> grad(const Point<dim> &point) const final

Gradient evaluator operator.

Parameters:

point – Point at which the function’s gradient is evaluated.

Returns:

Function gradient at point.

virtual Gradient<Interval<dim>> grad(const Point<dim, Interval<dim>> &point) const final

Gradient evaluator operator.

Parameters:

point – Point at which the function’s gradient is evaluated.

Returns:

Function gradient at point.

virtual Hessian<real> hessian(const Point<dim> &point) const final

Hessian evaluator operator.

Parameters:

point – Point at which the function’s hessian is evaluated.

Returns:

Function hessian at point.

template<int dim_aux = dim>
std::shared_ptr<MonomialsTP<dim - 1, range>> extract_facet(const int local_facet_id) const

Extracts a facet from the monomials tensor product.

This function extracts a facet of the monomials tensor product specified by the local facet ID.

Template Parameters:
  • dim – The dimension of the monomials tensor product.

  • range – The range type of the monomials tensor product.

Parameters:

local_facet_id – The ID of the local facet to extract.

Returns:

A shared pointer to the extracted facet.

Public Static Functions

template<typename T>
static Value<T> horner(const Point<dim, T> &point, typename std::vector<CoefsType>::const_iterator &coefs, const Vector<int, dim> &order)

Evaluates a monomials using Horner’s method.

For the method details check: https://en.wikipedia.org/wiki/Horner%27s_method

The evaluation along the dimensions is performed by recursively calling this function.

Template Parameters:

T – Type of the input and output variables.

Parameters:
  • point – Evaluation point.

  • coefs – Iterator to the coefficients of the polynomial. Their ordering is the same as the coefficients members of the class.

  • order – Order of the polynomial along the dim dimensions.

Returns:

Value of the polynomial at point.

template<typename T>
static Vector<Value<T>, dim + 1> horner_der(const Point<dim, T> &point, typename std::vector<CoefsType>::const_iterator &coefs, const Vector<int, dim> &order)

Evaluates the gradient of the monomials using the Horner’s method.

For the method details check: https://en.wikipedia.org/wiki/Horner%27s_method

The evaluation along the dimensions is performed by recursively calling this function.

Template Parameters:

T – Type of the input and output variables.

Parameters:
  • point – Evaluation point.

  • coefs – Iterator to the coefficients of the polynomial. Their ordering is the same as the coefficients members of the class.

  • order – Order of the polynomial along the dim dimensions.

Returns:

The first component of the return vector corresponds to the value of the polynomial itself, while the gradient is stored in the following dim components.

static void transform_coefs_to_Bezier(const MonomialsTP<dim, range> &monomials, std::vector<CoefsType> &bzr_coefs)

Transforms the coefficients of monomials to Bezier form.

This function takes the coefficients of monomials and transforms them into the corresponding Bezier coefficients.

Parameters:
  • monomials – Monomials whose coefficients are transformed.

  • bzr_coefs – A vector containing the Bezier coefficients to be computed.

template<int dim, int range>
class PolynomialTP : public qugar::impl::DomainFunc<dim, range>
#include <polynomial_tp.hpp>

Base class for tensor-product polynomial functions.

Template Parameters:
  • dim – Dimension of the parametric domain.

  • range – Dimension of the image.

Subclassed by qugar::impl::BezierTP< 2, 1 >, qugar::impl::BezierTP< dim, 1 >, qugar::impl::BezierTP< 3, 1 >, qugar::impl::BezierTP< dim, range >

Public Types

using CoefsType = std::conditional_t<range == 1, real, Point<range>>

Coefs type.

Public Functions

explicit PolynomialTP(const TensorSizeTP<dim> &order)

Constructor.

The coefficients vector is allocated by the constructor.

Parameters:

order – Order of the polynomial along each parametric direction.

PolynomialTP(const TensorSizeTP<dim> &order, const CoefsType &value)

Constructs a constant value PolynomialTP object with the specified order.

Parameters:
  • order – The order of the polynomial tensor.

  • value – Constant value for all the coefficients.

PolynomialTP(const TensorSizeTP<dim> &order, const std::vector<CoefsType> &coefs)

Constructor.

Parameters:
  • coefs – Coefficients of the polynomial. The ordering is such that dimension 0 is inner-most, i.e., iterates the fastest, while dimension dim-1 is outer-most and iterates the slowest.

  • order – Order of the polynomial along each parametric direction.

std::size_t get_num_coefs() const

Get the number of coeficients.

Returns:

Number of polynomial coefficients.

const std::vector<CoefsType> &get_coefs() const

Gets the polynomial coefficients.

Returns:

Constant reference to the polynomial coefficients.

const TensorSizeTP<dim> &get_order() const

Gets the polynomial order along the parametric directions.

Returns:

Constant reference to polynomial order.

int get_order(int dir) const

Gets the polynomial order (degree + 1) along the direction dir.

Parameters:

dir – Direction along which the order (degree + 1) is retrieved.

Returns:

Polynomial order (degree + 1) along dir.

const CoefsType &get_coef(int index) const

Retrieves the coefficients at the specified index.

This function returns a constant reference to the coefficients of the polynomial at the given index. The coefficients are of type CoefsType.

Parameters:

index – The index of the coefficients to retrieve.

Returns:

A constant reference to the coefficients at the specified index.

CoefsType &get_coef(int index)

Retrieves the coefficients at the specified index.

This function returns a constant reference to the coefficients of the polynomial at the given index. The coefficients are of type CoefsType.

Parameters:

index – The index of the coefficients to retrieve.

Returns:

A non-constant reference to the coefficients at the specified index.

const CoefsType &get_coef(const TensorIndexTP<dim> &index) const

Retrieves the coefficient associated with the given tensor index.

This function returns a reference to the coefficient corresponding to the specified tensor index. The coefficient is part of the polynomial tensor product representation.

Template Parameters:

dim – The dimension of the tensor index.

Parameters:

index – The tensor index for which the coefficient is to be retrieved.

Returns:

A constant reference to the coefficient associated with the given tensor index.

CoefsType &get_coef(const TensorIndexTP<dim> &index)

Retrieves the coefficient associated with the given tensor index.

This function returns a reference to the coefficient corresponding to the specified tensor index. The coefficient is part of the polynomial tensor product representation.

Template Parameters:

dim – The dimension of the tensor index.

Parameters:

index – The tensor index for which the coefficient is to be retrieved.

Returns:

A non-constant reference to the coefficient associated with the given tensor index.

int get_degree(int dir) const

Gets the polynomial degree (order - 1) along the direction dir.

Parameters:

dir – Direction along which the order is retrieved.

Returns:

Polynomial degree (order - 1) along dir.

void transform_image(const BoundBox<range> &old_domain, const BoundBox<range> &new_domain)

Transforms (in place) the coefficients of the polynomial from old_domain to new_domain.

Parameters:
  • old_domain – Domain from which polynomial coefficients are transformed from.

  • new_domain – Domain to which polynomial coefficients are transformed to.

void coefs_linear_transform(const real scale, const CoefsType &shift)

Applies a linear transformation to every coefficient.

It modifies every coefficient as coef = scale * coef + shift.

Parameters:
  • scale – Scale of every coefficient.

  • shift – Shift of every coefficient.

template<int dim>
class RefSystem
#include <ref_system.hpp>

A class representing a reference system in a given dimension.

This class provides functionality to create and manage a reference system in a specified dimension. It supports creating Cartesian reference systems centered at the origin with canonical directions, as well as custom reference systems with specified origins and axes.

Template Parameters:

dim – The dimension of the reference system.

Public Functions

RefSystem()

Constructs a new RefSystem object.

This is the default constructor for the RefSystem class.

It creates a reference system with the origin at the origin point and the basis with along the Cartesian axes.

explicit RefSystem(const Point<dim> &origin)

Constructs a RefSystem object with the specified origin point.

It creates a reference system with the origin at the given origin point and the basis with along the Cartesian axes.

Parameters:

origin – The origin point of the reference system.

RefSystem(const Point<dim> &origin, const Point<dim> &axis)

Constructs a reference system with a specified origin.

The given axis has different meanings depending on the dimension of the reference system. For 1D, the axis is the x-axis; for 2D, the axis is the x-axis; and for 3D, the axis is the z-axis. Thus, in 2D, the y-axis is computed by rotating the given axis 90 degrees. In 3D, given axis is considered to be the z-axis, then, x- and y-axes are computed to be orthonormal to axis.

Parameters:
  • origin – The origin point of the reference system.

  • axis – The axis point of the reference system.

template<int dim_aux = dim>
RefSystem(const Point<dim> &origin, const Point<dim> &axis_x, const Point<dim> &axis_y)

Constructs a reference system with a specified origin and two axes.

This constructor is only available for 3-dimensional reference systems.

The given axis_x and axis_y define the xy-plane, but they may not be orthonormal vectors. Thus, the way in which the new system directions are computed is as follows: First, the axis_x is normalized. Then, the z-axis is computed as the normalized cross-product of the axis_x and axis_y. Then, the y-axis is recomputed as the cross-product between normalized vectors along the z- and the x-axes.

Template Parameters:

dim_aux – Auxiliary template parameter to enforce the dimension check.

Parameters:
  • origin – The origin point of the reference system.

  • axis_x – The x-axis point of the reference system.

  • axis_y – The y-axis point of the reference system.

const Point<dim> &get_origin() const

Gets the origin point of the reference system.

Returns:

Constant reference to the origin point of the reference system.

const std::array<Point<dim>, dim> &get_basis() const

Gets the orthonormal basis of the reference system.

Returns:

Constant reference to the orthonormal basis of the reference system.

bool is_Cartesian_oriented() const

Checks if the reference system is Cartesian oriented.

This function determines whether the current reference system is oriented in a Cartesian manner. I.e., if the basis vectors are oriented along the Cartesian directions x, y, z.

Returns:

true if the reference system is Cartesian oriented, false otherwise.

template<int dim>
struct RootsIntervals
#include <impl_utils.hpp>

Struct for storing and managing computed roots and intervals of an implicit function.

Template Parameters:

dim – Dimension of the point at which the intervals are computed.

Public Functions

RootsIntervals()

Default constructor.

void clear()

Clears the container to the initial state.

void add_root(real root, int func_id)

Adds a new root.

Note

The roots are neither sorted nor adjusted (for degeneracies) after appending the new root.

Parameters:
  • root – New root to be added.

  • func_id – If of the restriction to which the root belongs to.

bool empty() const

Checks whether the container is empty.

Returns:

True if empty, i.e. there are no roots, false otherwise.

int get_num_roots() const

Gets the number of roots in the container.

Returns:

Number of roots.

void sort_roots()

Sorts (in increasing order) the roots in the container and the according restriction indices func_ids.

void adjust_roots(const Tolerance &tol, real x0, real x1)

Adjust the container roots by sorting them and forcing near roots (up to a _olerance) to be coincident.

Parameters:
  • tolTolerance to be used in the comparisons between roots.

  • x0 – Start of the interval to which the roots belong to.

  • x1 – End of the interval to which the roots belong to.

Public Members

std::vector<real> roots

List of roots.

Point<dim> point

Point at which roots are computed.

std::vector<int> func_ids

Restrictions to which root correspond to (there is a one to correspondence).

std::vector<bool> active_intervals

Flags indicating if the intervals defined by two consecutive roots are active.

template<int dim>
class UnfittedImplDomain : public qugar::UnfittedDomain<dim>
#include <impl_unfitted_domain.hpp>
namespace funcs

Namespace for defining implicit function examples. These function are ready to be consumed by Algoim.

template<int dim>
class AddFunctions : public qugar::impl::DomainFunc<dim>
#include <impl_funcs_lib.hpp>

This function adds two functions together.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit AddFunctions(const std::shared_ptr<const ImplicitFunc<dim>> &lhs, const std::shared_ptr<const ImplicitFunc<dim>> &rhs)

Constructor.

Parameters:
  • lhs – Left-hand-side operand.

  • rhs – Right-hand-side operand.

class Annulus : public qugar::impl::funcs::AnnulusBase, public qugar::impl::DomainFunc<2>
#include <primitive_funcs_lib.hpp>

2D annulus function. The function is defined by the annulus center and outer and inner radii. The function presents a negative sign inside the annulus (between both boundaries), and positive outside.

Note

non-Bezier version.

Public Functions

Annulus(real inner_radius, real outer_radius)

Constructor.

Center is set to (0,0).

Parameters:
  • inner_radius – Inner radius of the annulus.

  • outer_radius – Outer radius of the annulus.

Annulus(real inner_radius, real outer_radius, const Point<2> &center)

Constructor.

Parameters:
  • inner_radius – Inner radius of the annulus.

  • outer_radius – Outer radius of the annulus.

  • centerAnnulus’ center.

class AnnulusBase
#include <primitive_funcs_lib.hpp>

2D annulus base class. The function is defined by the annulus center and major and inner radii. The function presents a negative sign inside the annulus (between both boundaries), and positive outside.

Subclassed by qugar::impl::funcs::Annulus, qugar::impl::funcs::AnnulusBzr

Public Functions

AnnulusBase(real inner_radius, real outer_radius, const Point<2> &center)

Constructor.

Parameters:
  • inner_radius – Inner radius of the annulus.

  • outer_radius – Outer radius of the annulus.

  • centerAnnulus’ center.

const Point<2> &center() const

Gets the center of the annulus.

Returns:

real The annulus’s center.

real inner_radius() const

Gets the inner radius of the annulus.

Returns:

real The annulus’s inner radius.

real outer_radius() const

Gets the outer radius of the annulus.

Returns:

real The annulus’s outer radius.

class AnnulusBzr : public qugar::impl::funcs::AnnulusBase, public qugar::impl::BezierTP<2, 1>
#include <primitive_funcs_lib.hpp>

2D annulus function. The function is defined by the annulus center and outer and inner radii. The function presents a negative sign inside the annulus (between both boundaries), and positive outside.

Note

Bezier version.

Public Functions

AnnulusBzr(real inner_radius, real outer_radius)

Constructor.

Center is set to (0,0).

Parameters:
  • inner_radius – Inner radius of the annulus.

  • outer_radius – Outer radius of the annulus.

AnnulusBzr(real inner_radius, real outer_radius, const Point<2> &center)

Constructor.

Parameters:
  • inner_radius – Inner radius of the annulus.

  • outer_radius – Outer radius of the annulus.

  • centerAnnulus’ center.

template<int dim>
class Constant : public qugar::impl::funcs::ConstantBase, public qugar::impl::DomainFunc<dim>
#include <primitive_funcs_lib.hpp>

Dimension independent constant function.

Note

Non-Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

Constant()

Default constructor. Sets constant value to 0.5.

explicit Constant(real value)

Constructor.

Parameters:

valueConstant value.

class ConstantBase
#include <primitive_funcs_lib.hpp>

Dimension independent constant function.

Template Parameters:

dim – Parametric dimension.

Subclassed by qugar::impl::funcs::Constant< dim >, qugar::impl::funcs::ConstantBzr< dim >

Public Functions

explicit ConstantBase(real value)

Constructor.

Parameters:

valueConstant value.

real value() const

Gets the constant value.

Returns:

real The constant value.

template<int dim>
class ConstantBzr : public qugar::impl::funcs::ConstantBase, public qugar::impl::BezierTP<dim, 1>
#include <primitive_funcs_lib.hpp>

Dimension independent constant function.

Note

Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

ConstantBzr()

Default constructor. Sets constant value to 0.5.

explicit ConstantBzr(real value)

Constructor.

Parameters:

valueConstant value.

class Cylinder : public qugar::impl::funcs::CylinderBase, public qugar::impl::DomainFunc<3>
#include <primitive_funcs_lib.hpp>

Infinite cylinder.

The cylinder is defined by its radius, origin, and axis.

The function presents a negative sign around the cylinder’s axis, and positive far away. At a radius distance from the cylinder’s axis, the function vanishes.

Note

Non-Bezier version.

Public Functions

explicit Cylinder(real radius)

Constructor.

Cylinder along z-axis.

Parameters:

radiusCylinder’s radius.

Cylinder(real radius, const Point<3> &origin)

Constructor.

Cylinder with vertical axis at the given origin.

Parameters:
Cylinder(real radius, const Point<3> &origin, const Point<3> &axis)

Constructor.

Parameters:
class CylinderBase
#include <primitive_funcs_lib.hpp>

Infinite cylinder base class.

The cylinder is defined by its radius, origin, and axis.

The function presents a negative sign around the cylinder’s axis, and positive far away. At a radius distance from the cylinder’s axis, the function vanishes.

Subclassed by qugar::impl::funcs::Cylinder, qugar::impl::funcs::CylinderBzr

Public Functions

CylinderBase(real radius, const Point<3> &origin, const Point<3> &axis)

Constructor.

Parameters:
real radius() const

Gets the radius of the cylinder.

Returns:

real The cylinder’s radius.

const Point<3> &origin() const

Gets the origin of the cylinder.

Returns:

real The cylinder’s origin.

const Point<3> &axis() const

Gets the axis of the cylinder.

Returns:

real The cylinder’s axis.

class CylinderBzr : public qugar::impl::funcs::CylinderBase, public qugar::impl::BezierTP<3, 1>
#include <primitive_funcs_lib.hpp>

Infinite cylinder.

The cylinder is defined by its radius, origin, and axis.

The function presents a negative sign around the cylinder’s axis, and positive far away. At a radius distance from the cylinder’s axis, the function vanishes.

Note

Bezier version.

Public Functions

explicit CylinderBzr(real radius)

Constructor.

Cylinder along z-axis.

Parameters:

radiusCylinder’s radius.

CylinderBzr(real radius, const Point<3> &origin)

Constructor.

Cylinder with vertical axis at the given origin.

Parameters:
CylinderBzr(real radius, const Point<3> &origin, const Point<3> &axis)

Constructor.

Parameters:

Public Static Functions

static std::shared_ptr<MonomialsTP<3, 1>> create_monomials(real radius, const Point<3> &origin, const Point<3> &axis)

Creates a polynomial representation based on monomials for the given center, radius, and axis.

Parameters:
  • radius – The radius of the cylinder.

  • origin – The center of the cylinder.

  • axis – The radius of the cylinder.

template<int dim>
class DimLinear : public qugar::impl::funcs::FuncWithAffineTransf<dim>
#include <impl_funcs_lib.hpp>

dim-linear function.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit DimLinear(const std::array<real, num_coeffs> &coefs)

Constructs a new dim-linear function from its coefficients.

Parameters:

coefs – Function coefficients.

DimLinear(const std::array<real, num_coeffs> &coefs, const AffineTransf<dim> &transf)

Constructs a new dim-linear function from its coefficients.

Parameters:
  • coefs – Coefficients defining the function (stored in lexicographical ordering).

  • transf – Affine transformation applied to the dim-linear function. It may rotate the axes, translate the square, and/or scale it (iso or anisotropically).

Public Static Attributes

static const int num_coeffs = _impl::pow(2, dim)

Number of coefficients.

template<int dim>
class Ellipsoid : public qugar::impl::funcs::EllipsoidBase<dim>, public qugar::impl::DomainFunc<dim>
#include <primitive_funcs_lib.hpp>

Dimension independent ellipsoidal function. The function is defined by the ellipsoid’s semi-axes and centered at the origin. The function presents a negative sign around the origin, and positive far away.

Note

This class is implemented as an orthotropic scaling of a sphere.

Note

Non-Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit Ellipsoid(const Point<dim> &semi_axes)

Constructor.

Parameters:

semi_axes – Semi-axes length along the Cartesian axes.

Ellipsoid(const Point<dim> &semi_axes, const RefSystem<dim> &system)

Constructs an Ellipsoid object with specified semi-axes and reference system.

Parameters:
  • semi_axes – A Point object representing the semi-axes of the ellipsoid.

  • system – A RefSystem object representing the reference system in which the ellipsoid is defined.

template<int dim>
class EllipsoidBase
#include <primitive_funcs_lib.hpp>

Dimension independent ellipsoidal function (base cass). The function is defined by the ellipsoid’s semi-axes and centered at the origin. The function presents a negative sign around the origin, and positive far away.

Note

This class is implemented as an orthotropic scaling of a sphere.

Template Parameters:

dim – Parametric dimension.

Subclassed by qugar::impl::funcs::Ellipsoid< dim >, qugar::impl::funcs::EllipsoidBzr< dim >

Public Functions

EllipsoidBase(const Point<dim> &semi_axes, const RefSystem<dim> &system)

Constructs an Ellipsoid object with specified semi-axes and reference system.

Parameters:
  • semi_axes – A Point object representing the semi-axes of the ellipsoid.

  • system – A RefSystem object representing the reference system in which the ellipsoid is defined.

const Point<dim> &semi_axes() const

Gets the semi-axes of the ellipsoid.

Returns:

real The ellipsoid’s semi-axes.

const RefSystem<dim> &ref_system() const

Gets the reference system of the ellipsoid.

Returns:

real The ellipsoid’s reference system.

template<int dim>
class EllipsoidBzr : public qugar::impl::funcs::EllipsoidBase<dim>, public qugar::impl::BezierTP<dim, 1>
#include <primitive_funcs_lib.hpp>

Dimension independent ellipsoidal function. The function is defined by the ellipsoid’s semi-axes and centered at the origin. The function presents a negative sign around the origin, and positive far away.

Note

This class is implemented as an orthotropic scaling of a sphere.

Note

Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit EllipsoidBzr(const Point<dim> &semi_axes)

Constructor.

Parameters:

semi_axes – Semi-axes length along the Cartesian axes.

EllipsoidBzr(const Point<dim> &semi_axes, const RefSystem<dim> &system)

Constructs an Ellipsoid object with specified semi-axes and reference system.

Parameters:
  • semi_axes – A Point object representing the semi-axes of the ellipsoid.

  • system – A RefSystem object representing the reference system in which the ellipsoid is defined.

template<int dim>
class FuncWithAffineTransf : public qugar::impl::DomainFunc<dim>
#include <impl_funcs_lib.hpp>

Function containing an affine transformation.

Template Parameters:

dim – Parametric dimension.

Subclassed by qugar::impl::funcs::DimLinear< dim >, qugar::impl::funcs::Square< dim >, qugar::impl::funcs::TransformedFunction< dim >

Public Functions

FuncWithAffineTransf()

Default constructor. Creates and stores an identity transformation.

explicit FuncWithAffineTransf(const AffineTransf<dim> &transf)

Constructs a new class storing the given transf.

Parameters:

transf – Transformation to store.

template<int dim>
class Negative : public qugar::impl::DomainFunc<dim>
#include <impl_funcs_lib.hpp>

This function computes the negative of a given function.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit Negative(const std::shared_ptr<const ImplicitFunc<dim>> &func)

Constructor.

Parameters:

func – Function whose negative is computed.

template<int dim>
class Plane : public qugar::impl::funcs::PlaneBase<dim>, public qugar::impl::DomainFunc<dim>
#include <primitive_funcs_lib.hpp>

Plane function.

This is a linear function whose value is zero at a line, and grows linearly (positively or negatively) as you move away from the line.

Note

Non-Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

Plane()

Constructs a new plane function. The line (levelset) is the line x=0.

Plane(const Point<dim> &origin, const Point<dim> &normal)

Constructs a new plane function. The line is defined by an origin and a normal vector.

template<int dim>
class PlaneBase
#include <primitive_funcs_lib.hpp>

Plane base class.

This is a linear function whose value is zero at a line, and grows linearly (positively or negatively) as you move away from the line.

Template Parameters:

dim – Parametric dimension.

Subclassed by qugar::impl::funcs::Plane< dim >, qugar::impl::funcs::PlaneBzr< dim >

Public Functions

PlaneBase(const Point<dim> &origin, const Point<dim> &normal)

Constructs a new plane function. The line is defined by an origin and a normal vector.

const Point<dim> &origin() const

Gets the origin of the plane.

Returns:

real The plane’s origin.

const Point<dim> &normal() const

Gets the normal of the plane.

Returns:

real The plane’s normal.

template<int dim>
class PlaneBzr : public qugar::impl::funcs::PlaneBase<dim>, public qugar::impl::BezierTP<dim, 1>
#include <primitive_funcs_lib.hpp>

Plane function.

This is a linear function whose value is zero at a line, and grows linearly (positively or negatively) as you move away from the line.

Note

Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

PlaneBzr()

Constructs a new plane function. The line (levelset) is the line x=0.

PlaneBzr(const Point<dim> &origin, const Point<dim> &normal)

Constructs a new plane function. The line is defined by an origin and a normal vector.

template<int dim>
class Sphere : public qugar::impl::funcs::SphereBase<dim>, public qugar::impl::DomainFunc<dim>
#include <primitive_funcs_lib.hpp>

Dimension independent spherical function. The function is defined by its center and radius. The function presents a negative sign around the center, and positive far away. At a radius distance from the center, the function vanishes.

Note

Non-Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit Sphere(real radius)

Constructs a Sphere with the given radius and centered at the origin.

Parameters:

radius – The radius of the sphere.

Sphere(real radius, const Point<dim> &center)

Constructs a Sphere object with a specified center and radius.

Parameters:
  • radius – The radius of the sphere.

  • center – The center of the sphere.

template<int dim>
class SphereBase
#include <primitive_funcs_lib.hpp>

Dimension independent spherical function.

Note

This is the base class for Bezier and general implementations.

Template Parameters:

dim – Parametric dimension.

Subclassed by qugar::impl::funcs::Sphere< dim >, qugar::impl::funcs::SphereBzr< dim >

Public Functions

SphereBase(real radius, const Point<dim> &center)

Constructs a Sphere object with a specified center and radius.

Parameters:
  • radius – The radius of the sphere.

  • center – The center of the sphere.

real radius() const

Gets the radius of the sphere.

Returns:

real The sphere’s radius.

const Point<dim> &center() const

Gets the center of the sphere.

Returns:

real The sphere’s center.

template<int dim>
class SphereBzr : public qugar::impl::funcs::SphereBase<dim>, public qugar::impl::BezierTP<dim, 1>
#include <primitive_funcs_lib.hpp>

Dimension independent spherical function. The function is defined by its center and radius. The function presents a negative sign around the center, and positive far away. At a radius distance from the center, the function vanishes.

Note

Bezier version.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit SphereBzr(real radius)

Constructs a Sphere with the given radius and centered at the origin.

Parameters:

radius – The radius of the sphere.

SphereBzr(real radius, const Point<dim> &center)

Constructs a Sphere object with a specified center and radius.

Parameters:
  • radius – The radius of the sphere.

  • center – The center of the sphere.

template<int dim>
class Square : public qugar::impl::funcs::FuncWithAffineTransf<dim>
#include <impl_funcs_lib.hpp>

Dimension independent square function in domain [-1,1]^dim.

Warning

So far, only implemented for the 2D case.

Template Parameters:

dim – Parametric dimension.

Public Functions

Square() = default

Constructs a new square aligned with the Cartesian axes.

explicit Square(const AffineTransf<dim> &transf)

Constructs a new square object transformed according to transf.

Parameters:

transf – Affine transformation applied to the square. It may rotate the axes, translate the square, and/or scale it (iso or anisotropically).

Public Static Functions

template<typename T>
static int sgn(const T &val)

Returns the sign of the given value.

Template Parameters:

T – Type of the value.

Parameters:

val – Value to be tested.

Returns:

+1 is val is positive, -1 if it is negative, 0 otherwise.

template<int dim>
class SubtractFunctions : public qugar::impl::DomainFunc<dim>
#include <impl_funcs_lib.hpp>

This function subtracts two functions.

Template Parameters:

dim – Parametric dimension.

Public Functions

explicit SubtractFunctions(const std::shared_ptr<const ImplicitFunc<dim>> &lhs, const std::shared_ptr<const ImplicitFunc<dim>> &rhs)

Constructor.

Parameters:
  • lhs – Left-hand-side operand.

  • rhs – Right-hand-side operand.

class Torus : public qugar::impl::funcs::TorusBase, public qugar::impl::DomainFunc<3>
#include <primitive_funcs_lib.hpp>

3D torus function. The function is defined by the torus center, axis, and major and inner radii. The function presents a negative sign inside the torus, and positive outside.

Note

Non-Bezier version.

Public Functions

Torus(real major_radius, real minor_radius)

Constructor.

Center is set to (0,0,0) and axis to (0,0,1).

Parameters:
  • major_radius – Radius of the major circle.

  • minor_radius – Radius of the inner circle.

Torus(real major_radius, real minor_radius, const Point<3> &center)

Constructor.

The axis is set to (0,0,1).

Parameters:
  • major_radius – Radius of the major circle.

  • minor_radius – Radius of the inner circle.

  • centerTorus’ center.

Torus(real major_radius, real minor_radius, const Point<3> &center, const Point<3> &axis)

Constructor.

Parameters:
  • major_radius – Radius of the major circle.

  • minor_radius – Radius of the inner circle.

  • centerTorus’ center.

  • axisTorus’ axis (perpendicular to the major circle plane).

template<typename T>
Point<3, T> compute_P_x_0(const Point<3, T> &point) const

Computes the normal component of the a vector respect to the origin.

Given a point, computes the relative vector respect to the origin, and then computes the normal component of the vector respect to the torus’ axis.

Template Parameters:

T – The type of the coordinates of the point.

Parameters:

point – The input 3D point for which the P_x_0 value is to be computed.

Returns:

A 3D point representing the computed P_x_0 value.

class TorusBase
#include <primitive_funcs_lib.hpp>

3D torus function base class. The function is defined by the torus center, axis, and major and inner radii. The function presents a negative sign inside the torus, and positive outside.

Subclassed by qugar::impl::funcs::Torus, qugar::impl::funcs::TorusBzr

Public Functions

TorusBase(real major_radius, real minor_radius, const Point<3> &center, const Point<3> &axis)

Constructor.

Parameters:
  • major_radius – Radius of the major circle.

  • minor_radius – Radius of the inner circle.

  • centerTorus’ center.

  • axisTorus’ axis (perpendicular to the major circle plane).

real major_radius() const

Gets the major radius of the torus.

Returns:

real The torus’ major radius.

real minor_radius() const

Gets the minor radius of the torus.

Returns:

real The torus’ minor radius.

const Point<3> &center() const

Gets the center of the plane.

Returns:

real The plane’s center.

const Point<3> &axis() const

Gets the axis of the plane.

Returns:

real The plane’s axis.

class TorusBzr : public qugar::impl::funcs::TorusBase, public qugar::impl::BezierTP<3, 1>
#include <primitive_funcs_lib.hpp>

3D torus function. The function is defined by the torus center, axis, and major and inner radii. The function presents a negative sign inside the torus, and positive outside.

Note

Bezier version.

Public Functions

TorusBzr(real major_radius, real minor_radius)

Constructor.

Center is set to (0,0,0) and axis to (0,0,1).

Parameters:
  • major_radius – Radius of the major circle.

  • minor_radius – Radius of the inner circle.

TorusBzr(real major_radius, real minor_radius, const Point<3> &center)

Constructor.

The axis is set to (0,0,1).

Parameters:
  • major_radius – Radius of the major circle.

  • minor_radius – Radius of the inner circle.

  • centerTorus’ center.

TorusBzr(real major_radius, real minor_radius, const Point<3> &center, const Point<3> &axis)

Constructor.

Parameters:
  • major_radius – Radius of the major circle.

  • minor_radius – Radius of the inner circle.

  • centerTorus’ center.

  • axisTorus’ axis (perpendicular to the major circle plane).

template<int dim>
class TransformedFunction : public qugar::impl::funcs::FuncWithAffineTransf<dim>
#include <impl_funcs_lib.hpp>

Creates a new implicit function that is just a base function to which an affine transformation is applied.

Template Parameters:

dim – Function’s parametric direction.

Public Functions

TransformedFunction(const std::shared_ptr<const ImplicitFunc<dim>> &base_func, const AffineTransf<dim> &transf)

Constructs a new transformed function.

Parameters:
  • base_func – Base function to be transformed.

  • transf – Affine transformation to apply.

namespace _impl

Typedefs

template<typename V, int new_dim>
using NewVector_t = typename NewVector<V, new_dim>::type
template<typename V>
using Hessian_t = typename NewVector<V, (VectorDim<V>::dim * (VectorDim<V>::dim + 1)) / 2>::type
template<typename V>
using VectorType_t = typename VectorType<V>::type

Functions

template<class T>
T pow(const T base, unsigned const exponent) noexcept

Computes the exponent of a value.

Template Parameters:

T – Input and output’s type.

Parameters:
  • base – Base of the exponent.

  • exponent – Exponent of the base.

Returns:

Computed value

template<typename V>
struct IsAlgoimVector
template<int dim>
struct IsAlgoimVector<Point<dim>>
#include <impl_funcs_lib.hpp>
template<typename T, int dim>
struct IsAlgoimVector<Vector<T, dim>>
#include <impl_funcs_lib.hpp>
template<typename V, int new_dim>
struct NewVector
template<int dim, int new_dim>
struct NewVector<Point<dim>, new_dim>
#include <impl_funcs_lib.hpp>
template<typename T, int dim, int new_dim>
struct NewVector<Vector<T, dim>, new_dim>
#include <impl_funcs_lib.hpp>
template<typename V>
struct VectorDim
template<int dim_>
struct VectorDim<Point<dim_>>
#include <impl_funcs_lib.hpp>
template<typename T, int dim_>
struct VectorDim<Vector<T, dim_>>
#include <impl_funcs_lib.hpp>
template<typename V>
struct VectorType
template<int dim>
struct VectorType<Point<dim>>
#include <impl_funcs_lib.hpp>
template<typename T, int dim>
struct VectorType<Vector<T, dim>>
#include <impl_funcs_lib.hpp>
namespace tpms

Namespace for Triple-Periodic Minimal Surfaces. Namely, Schoen gyroid, Schoen IWP, Scheon FRD, Fischer-Koch S, Schwarz diamond, and Schwarz primitive. These function are ready to be consumed by Algoim.

Functions

declare_tpms(Schoen)

Schoen’s gyroid function. Defined as f(x,y,z,m,n,q) = sin(2 pi m x) * cos(2 pi n y) + sin(2 pi n y) * cos(2 pi q z) + sin(2 pi q z) * cos(2 pi m x) this is a triply periodic function with period (m, n, q).

See https:// en.wikipedia.org/wiki/Gyroid

declare_tpms(SchoenIWP)

Schoen IWP’s gyroid function. Defined as f(x,y,z,m,n,q) = 2 * (cos(2 pi m x) * cos(2 pi n y) + cos(2 pi n y) * cos(2 pi q z)

  • cos(2 pi q z) * cos(2 pi m x)) - cos(4 pi m x) - cos(4 pi m n y) - cos(4 pi q z) this is a triply periodic function with period (m, n, q).

See https://en.wikipedia.org/wiki/Gyroid

declare_tpms(SchoenFRD)

Schoen FRD’s gyroid function. Defined as f(x,y,z,m,n,q) = 4 * cos(2 pi m x) * cos(2 pi n y) * cos(2 pi q z) - cos(4 pi m x) * cos(4 pi n y)

  • cos(4 pi n y) * cos(4 pi q z) - cos(4 pi m q z* cos(4 pi m x) this is a triply periodic function with period (m, n, q).

See https://en.wikipedia.org/wiki/Gyroid

declare_tpms(FischerKochS)

Fischer-Koch S’ gyroid function. Defined as f(x,y,z,m,n,q) = cos(4 pi m x) * sin(2 pi n y) * cos(2 pi q z)

  • cos(2 pi m x) * cos(4 pi n y) * sin(2 pi q z)

  • sin(2 pi m x) * cos(2 pi n y) * cos(4 pi q z) this is a triply periodic function with period (m, n, q).

See https://en.wikipedia.org/wiki/Gyroid

declare_tpms(SchwarzDiamond)

Schwarz diamond’s gyroid function. Defined as f(x,y,z,m,n,q) = cos(2 pi m x) * cos(2 pi n y) * cos(2 pi q z)

  • sin(2 pi m x) * sin(2 pi n y) * sin(2 pi q z) this is a triply periodic function with period (m, n, q).

See https://en.wikipedia.org/wiki/Gyroid

declare_tpms(SchwarzPrimitive)

Schwarz primitive’s gyroid function. Defined as f(x,y,z,m,n,q) = cos(2 pi mpoint(0)) + cos(2 pi n y) + cos(2 pi q z) this is a triply periodic function with period (m, n, q).

See https://en.wikipedia.org/wiki/Gyroid

template<int dim>
class TPMSBase : public qugar::impl::DomainFunc<dim>
#include <tpms_lib.hpp>
namespace numbers

Numbers’ namespace.

Variables

real zero = {0.0}

Real zero value.

real one = {1.0}

Real one value.

real two = {2.0}

Real two value.

real three = {3.0}

Real three value.

real four = {4.0}

Real four value.

real five = {5.0}

Real five value.

real six = {6.0}

Real six value.

real seven = {7.0}

Real seven value.

real eight = {8.0}

Real eight value.

real nine = {9.0}

Real nine value.

real half = {0.5}

Real one over two value.

real one_third = {1.0 / 3.0}

Real one third value.

real two_thirds = {2.0 / 3.0}

Real two thirds value.

real four_thirds = {4.0 / 3.0}

Real four thirds value.

real one_quarter = {0.25}

Real one quarter value.

real three_quarters = {0.75}

Real three quarters value.

real pi = {std::numbers::pi_v<real>}

Pi.

real infty = {std::numeric_limits<real>::infinity()}

Infinity.

real eps = std::numeric_limits<real>::epsilon()

Machine epsilon.

real near_eps = real{10.0} * eps

Near machine epsilon (10 times the machine precision.