IgANet
IgANets - Isogeometric Analysis Networks
Loading...
Searching...
No Matches
iganet::UniformBSplineCore< real_t, GeoDim, Degrees > Class Template Reference

Tensor-product uniform B-spline (core functionality) More...

#include </home/runner/work/iganet/iganet/include/bspline.hpp>

Inheritance diagram for iganet::UniformBSplineCore< real_t, GeoDim, Degrees >:
iganet::utils::Serializable iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>

Public Types

template<typename other_t , short_t GeoDim_, short_t... Degrees_>
using derived_self_type = UniformBSplineCore< other_t, GeoDim_, Degrees_... >
 Deduces the derived self-type when exposed to different class template parameters real_t and GeoDim, and the Degrees parameter pack.
 
template<template< typename, short_t, short_t... > class BSpline, std::make_signed< short_t >::type degree_elevate = 0>
using derived_type = BSpline< real_t, GeoDim,(Degrees+degree_elevate)... >
 Deduces the type of the template template parameter BSpline when exposed to the class template parameters real_t and GeoDim, and the Degrees parameter pack. The optional template parameter degree_elevate can be used to (de-)elevate the degrees by an additive constant.
 
template<typename other_t >
using real_derived_self_type = UniformBSplineCore< other_t, GeoDim, Degrees... >
 Deduces the derived self-type when exposed to a different class template parameter real_t
 
template<std::make_signed< short_t >::type degree_elevate = 0>
using self_type = derived_type< UniformBSplineCore, degree_elevate >
 Deduces the self-type possibly degrees (de-)elevated by the additive constant degree_elevate
 
using value_type = real_t
 Value type.
 

Public Member Functions

 UniformBSplineCore (const std::array< int64_t, parDim_ > &ncoeffs, const utils::TensorArray< geoDim_ > &coeffs, bool clone=false, Options< real_t > options=Options< real_t >{})
 Constructor for equidistant knot vectors.
 
 UniformBSplineCore (const std::array< int64_t, parDim_ > &ncoeffs, enum init init=init::greville, Options< real_t > options=Options< real_t >{})
 Constructor for equidistant knot vectors.
 
 UniformBSplineCore (const std::array< int64_t, parDim_ > &ncoeffs, utils::TensorArray< geoDim_ > &&coeffs, Options< real_t > options=Options< real_t >{})
 Constructor for equidistant knot vectors.
 
template<typename other_t >
 UniformBSplineCore (const UniformBSplineCore< other_t, GeoDim, Degrees... > &other, Options< real_t > options=Options< real_t >{})
 Copy constructor.
 
 UniformBSplineCore (Options< real_t > options=Options< real_t >{})
 Default constructor.
 
torch::Tensor as_tensor () const noexcept override
 Returns all coefficients as a single tensor.
 
int64_t as_tensor_size () const noexcept override
 Returns the size of the single tensor representation of all coefficients.
 
const utils::TensorArray< geoDim_ > & coeffs () const noexcept
 Returns a constant reference to the array of coefficient vectors.
 
utils::TensorArray< geoDim_ > & coeffs () noexcept
 Returns a non-constant reference to the array of coefficient vectors.
 
const torch::Tensor & coeffs (short_t i) const noexcept
 Returns a constant reference to the coefficient vector in the \(i\)-th dimension.
 
torch::Tensor & coeffs (short_t i) noexcept
 Returns a non-constant reference to the coefficient vector in the \(i\)-th dimension.
 
nlohmann::json coeffs_to_json () const
 Returns the B-spline object's coefficients as JSON object.
 
utils::TensorArray< geoDim_coeffs_view () const noexcept
 Returns an array of views to the coefficient vectors.
 
const auto coeffs_view (short_t i) const noexcept
 Returns a view to the coefficient vector in the \(i\)-th dimension.
 
torch::Device device () const noexcept override
 Returns the device property.
 
int32_t device_index () const noexcept override
 Returns the device_index property.
 
torch::Dtype dtype () const noexcept override
 Returns the dtype property.
 
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto eval (const utils::TensorArray< parDim_ > &xi, const utils::TensorArray< parDim_ > &knot_indices) const
 Returns the value of the univariate B-spline object in the points xi
 
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto eval (const utils::TensorArray< parDim_ > &xi, const utils::TensorArray< parDim_ > &knot_indices, const torch::Tensor &coeff_indices) const
 Returns the value of the univariate B-spline object in the points xi
 
template<typename BSpline >
autofrom_gismo (BSpline &bspline, bool, bool)
 
UniformBSplineCorefrom_json (const nlohmann::json &json)
 Updates the B-spline object from JSON object.
 
UniformBSplineCorefrom_tensor (const torch::Tensor &tensor) noexcept override
 Sets all coefficients from a single tensor.
 
UniformBSplineCorefrom_xml (const pugi::xml_document &doc, int id=0, std::string label="", int index=-1)
 Updates the B-spline object from XML object.
 
UniformBSplineCorefrom_xml (const pugi::xml_node &root, int id=0, std::string label="", int index=-1)
 Updates the B-spline object from XML node.
 
auto greville (bool interior=false) const
 Returns the Greville abscissae.
 
void init_coeffs (enum init init)
 Initializes the B-spline coefficients.
 
void init_knots ()
 Initializes the B-spline knots.
 
bool is_sparse () const noexcept override
 Returns true if the layout is sparse.
 
template<typename other_t , short_t GeoDim_, short_t... Degrees_>
bool isclose (const UniformBSplineCore< other_t, GeoDim_, Degrees_... > &other, real_t rtol=real_t{1e-5}, real_t atol=real_t{1e-8}) const
 Returns true if both B-spline objects are close up to the given tolerances.
 
const utils::TensorArray< parDim_ > & knots () const noexcept
 Returns a constant reference to the array of knot vectors.
 
utils::TensorArray< parDim_ > & knots () noexcept
 Returns a non-constant reference to the array of knot vectors.
 
const torch::Tensor & knots (short_t i) const noexcept
 Returns a constant reference to the knot vector in the \(i\)-th dimension.
 
torch::Tensor & knots (short_t i) noexcept
 Returns a non-constant reference to the knot vector in the \(i\)-th dimension.
 
nlohmann::json knots_to_json () const
 Returns the B-spline object's knots as JSON object.
 
torch::Layout layout () const noexcept override
 Returns the layout property.
 
void load (const std::string &filename, const std::string &key="bspline")
 Loads the B-spline from file.
 
const std::array< int64_t, parDim_ > & ncoeffs () const noexcept
 Returns a constant reference to the array of coefficient vector dimensions.
 
int64_t ncoeffs (short_t i) const noexcept
 Returns the total number of coefficients in the \(i\)-th direction.
 
int64_t ncumcoeffs () const noexcept
 Returns the total number of coefficients.
 
const std::array< int64_t, parDim_ > & nknots () const noexcept
 Returns a constant reference to the array of knot vector dimensions.
 
int64_t nknots (short_t i) const noexcept
 Returns the dimension of the knot vector in the \(i\)-th dimension.
 
template<typename other_t , short_t GeoDim_, short_t... Degrees_>
bool operator!= (const UniformBSplineCore< other_t, GeoDim_, Degrees_... > &other) const
 Returns true if both B-spline objects are different.
 
template<typename other_t , short_t GeoDim_, short_t... Degrees_>
bool operator== (const UniformBSplineCore< other_t, GeoDim_, Degrees_... > &other) const
 Returns true if both B-spline objects are the same.
 
const Options< real_t > & options () const noexcept
 Returns a constant reference to the B-spline object's options.
 
bool pinned_memory () const noexcept override
 Returns the pinned_memory property.
 
torch::serialize::InputArchive & read (torch::serialize::InputArchive &archive, const std::string &key="bspline")
 Reads the B-spline from a torch::serialize::InputArchive object.
 
bool requires_grad () const noexcept override
 Returns the requires_grad property.
 
void save (const std::string &filename, const std::string &key="bspline") const
 Saves the B-spline to file.
 
UniformBSplineCoreset_requires_grad (bool requires_grad) noexcept override
 Sets the B-spline object's requires_grad property.
 
auto to_gismo () const
 Converts the B-spline object into a gsBSpline object of the parametric dimension is one and a gsTensorBSpline object otherwise.
 
template<typename BSpline >
BSplineto_gismo (BSpline &bspline, bool, bool) const
 
nlohmann::json to_json () const override
 Returns the B-spline object as JSON object.
 
pugi::xml_document to_xml (int id=0, std::string label="", int index=-1) const
 Returns the B-spline object as XML object.
 
pugi::xml_node & to_xml (pugi::xml_node &root, int id=0, std::string label="", int index=-1) const
 Returns the B-spline object as XML node.
 
UniformBSplineCoretransform (const std::function< std::array< real_t, geoDim_ >(const std::array< real_t, parDim_ > &)> transformation)
 Transforms the coefficients based on the given mapping.
 
template<std::size_t N>
UniformBSplineCoretransform (const std::function< std::array< real_t, N >(const std::array< real_t, parDim_ > &)> transformation, std::array< short_t, N > dims)
 Transforms the coefficients based on the given mapping.
 
UniformBSplineCoreuniform_refine (int numRefine=1, int dim=-1)
 Returns the B-spline object with uniformly refined knot and coefficient vectors.
 
torch::serialize::OutputArchive & write (torch::serialize::OutputArchive &archive, const std::string &key="bspline") const
 Writes the B-spline into a torch::serialize::OutputArchive object.
 
- Public Member Functions inherited from iganet::utils::Serializable
virtual void pretty_print (std::ostream &os=Log(log::info)) const =0
 Returns a string representation of the object.
 
- Public Member Functions inherited from iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>
virtual void pretty_print (std::ostream &os=Log(log::info)) const noexcept=0
 Returns a string representation.
 

Static Public Member Functions

static constexpr const short_tdegree (short_t i) noexcept
 Returns a constant reference to the degree in the \(i\)-th dimension.
 
static constexpr const std::array< short_t, parDim_ > & degrees () noexcept
 Returns a constant reference to the array of degrees.
 
static constexpr short_t geoDim () noexcept
 Returns the geometric dimension.
 
static constexpr bool is_nonuniform () noexcept
 Returns true if the B-spline is non-uniform.
 
static constexpr bool is_uniform () noexcept
 Returns true if the B-spline is uniform.
 
static constexpr short_t parDim () noexcept
 Returns the parametric dimension.
 

Protected Member Functions

template<short_t degree, short_t dim, short_t deriv>
auto eval_basfunc_univariate (const torch::Tensor &xi, const torch::Tensor &knot_indices) const
 Returns the vector of univariate B-spline basis functions (or their derivatives) evaluated in the point xi
 
void update_coeffs (const utils::TensorArray< parDim_ > &knots, const utils::TensorArray< parDim_ > &knot_indices)
 Updates the B-spline coefficients after knot insertion.
 
template<short_t degree, short_t dim>
auto update_coeffs_univariate (const torch::Tensor &knots, const torch::Tensor &knot_indices) const
 Returns the knot insertion matrix.
 

Protected Attributes

utils::TensorArray< geoDim_coeffs_
 Array storing the coefficients of the control net \(\left(\mathbf{c}_{i_d}\right)_{i_d=1}^{n_d}\), \(\mathbf{c}_{i_d}\in\mathbb{R}^{d_\text{geo}}\).
 
utils::TensorArray< parDim_knots_
 Array storing the knot vectors \(\left(\left(t_{i_d}\right)_{i_d=1}^{n_d+p_d+1}\right)_{d=1}^{d_\text{par}}\).
 
std::array< int64_t, parDim_ncoeffs_
 Array storing the sizes of the coefficients of the control net \(\left(n_d\right)_{d=1}^{d_\text{par}}\).
 
std::array< int64_t, parDim_ncoeffs_reverse_
 Array storing the sizes of the coefficients of the control net \(\left(n_d\right)_{d=1}^{d_\text{par}}\) in reverse order (needed for coeffs_view)
 
std::array< int64_t, parDim_nknots_
 Array storing the sizes of the knot vectors \(\left(n_d+p_d+1\right)_{d=1}^{d_\text{par}}\).
 
Options< real_t > options_
 Options.
 

Static Protected Attributes

static constexpr const std::array< short_t, parDim_degrees_ = {Degrees...}
 Array storing the degrees \(\left(p_d\right)_{d=1}^{d_\text{par}}\).
 
static constexpr const short_t geoDim_ = GeoDim
 Dimension of the geometric space \(\Omega\subset\mathbb{R}^{d_\text{geo}}\).
 
static constexpr const short_t parDim_ = sizeof...(Degrees)
 Dimension of the parametric space \(\hat\Omega=[0,1]^{d_\text{par}}\).
 

Private Member Functions

template<std::size_t... Is>
torch::Tensor as_tensor_ (std::index_sequence< Is... >) const noexcept
 Returns all coefficients as a single tensor.
 
template<int64_t degree, int64_t deriv, int64_t terminal = degree - deriv>
int64_t constexpr eval_prefactor () const
 Computes the prefactor \(p_d!/(p_d-r_d)! = p_d \cdots (p_d-r_d+1)\).
 
template<std::size_t... Is>
UniformBSplineCorefrom_tensor_ (std::index_sequence< Is... >, const torch::Tensor &tensor) noexcept
 Sets all coefficients from a single tensor.
 

Friends

template<typename BSplineCore >
class BSplineCommon
 Enable access to private members.
 

Detailed Description

template<typename real_t, short_t GeoDim, short_t... Degrees>
class iganet::UniformBSplineCore< real_t, GeoDim, Degrees >

Tensor-product uniform B-spline (core functionality)

This class implements the core functionality of all B-spline classes and serves as base class for (non-)uniform B-splines.

Mathematically, this class defines a mapping

\[ \mathbf{f}:\hat\Omega \mapsto \Omega \]

from the \(d_\text{par}\)-dimensional parametric space \(\hat\Omega=[0,1]^{d_\text{par}}\) to the \(d_\text{geo}\)-dimensional geometric space \(\Omega\subset\mathbb{R}^{d_\text{geo}}\).

This mapping is defined by tensor-product B-spline basis functions

\[ B_I(\boldsymbol{\xi}) = \bigotimes_{d=1}^{d_\text{par}} B_{i_d,p_d}(\xi_d) \]

and the control points

\[ \mathbf{c}_I = \mathbf{c}_{i_1,i_2,\dots, i_{d_\text{par}}} \in \mathbb{R}^{d_\text{geo}}. \]

Here, \(i_d\) are the local numbers of the univariate B-splines \(\left(B_{i_d,p_d}\right)_{i_d=1}^{n_d}\) in the \(d\)-th parametric dimension, \(p_d\) is the respective degree, and \(n_d\) is the number of univariate B-splines in the \(d\)-th direction. Moreover, \(0\le \xi_{i_d}\le 1\) is the parametric value at which the B-spline is evaluated. The multivariate B-spline function is defined as follows

\[ \mathbf{f}(\boldsymbol{\xi}) = \sum_{I=1}^N B_I(\boldsymbol{\xi}) \mathbf{c}_I \]

Here and below we adopt the vector notation \(\boldsymbol{\xi} = \left(\xi_1,\xi_2,\dots,\xi_{d_\text{par}}\right)^\top\) and combine multiple local indices \(i_1,i_2,\dots,i_{d_\text{par}}\) of univariate B-spline basis functions into the global index \(1\le I \le N\) with \(N=n_1\cdot n_2\cdot\dots\cdot n_{d_\text{par}}\) denoting the total number of multivariate B-splines.

This class implements B-spline functions and their derivatives for 1, 2, 3, and 4 parametric dimensions. The univariate B-splines are uniquely determined by their knot vectors

\[ \left(t_{i_d}\right)_{i_d=1}^{n_d+p_d+1} \]

with \(0\le t_{i_d}\le 1\) and \(t_{i_d}\le t_{i_d+1}\) for all \(i_d\), that is, the knot vectors are given by a non-decreasing sequence of values in the interval \([0,1]\) with the possibility that knot values are repeated.

This class implements the evaluation of B-splines and their derivatives as explained in Chapters 2 and 3 from [2].

Note
C++ uses 0-based indexing so that all of the above formulas need to be shifted by -1. Moreover, all vectors, matrices, and tensors are implemented as torch::Tensor objects and hence adopt Torch's local-to-global mapping. It is therefore imperative to always use Torch's indexing functionality to extract sub-tensors.

Member Typedef Documentation

◆ derived_self_type

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename other_t , short_t GeoDim_, short_t... Degrees_>
using iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::derived_self_type = UniformBSplineCore<other_t, GeoDim_, Degrees_...>

Deduces the derived self-type when exposed to different class template parameters real_t and GeoDim, and the Degrees parameter pack.

◆ derived_type

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<template< typename, short_t, short_t... > class BSpline, std::make_signed< short_t >::type degree_elevate = 0>
using iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::derived_type = BSpline<real_t, GeoDim, (Degrees + degree_elevate)...>

Deduces the type of the template template parameter BSpline when exposed to the class template parameters real_t and GeoDim, and the Degrees parameter pack. The optional template parameter degree_elevate can be used to (de-)elevate the degrees by an additive constant.

◆ real_derived_self_type

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename other_t >
using iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::real_derived_self_type = UniformBSplineCore<other_t, GeoDim, Degrees...>

Deduces the derived self-type when exposed to a different class template parameter real_t

◆ self_type

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<std::make_signed< short_t >::type degree_elevate = 0>
using iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::self_type = derived_type<UniformBSplineCore, degree_elevate>

Deduces the self-type possibly degrees (de-)elevated by the additive constant degree_elevate

◆ value_type

template<typename real_t , short_t GeoDim, short_t... Degrees>
using iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::value_type = real_t

Value type.

Constructor & Destructor Documentation

◆ UniformBSplineCore() [1/5]

template<typename real_t , short_t GeoDim, short_t... Degrees>
iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::UniformBSplineCore ( Options< real_t >  options = Options<real_t>{})
inline

Default constructor.

Parameters
[in]optionsOptions configuration

◆ UniformBSplineCore() [2/5]

template<typename real_t , short_t GeoDim, short_t... Degrees>
iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::UniformBSplineCore ( const std::array< int64_t, parDim_ > &  ncoeffs,
enum init  init = init::greville,
Options< real_t >  options = Options<real_t>{} 
)
inline

Constructor for equidistant knot vectors.

Parameters
[in]ncoeffsNumber of coefficients per parametric dimension
[in]initType of initialization
[in]optionsOptions configuration

◆ UniformBSplineCore() [3/5]

template<typename real_t , short_t GeoDim, short_t... Degrees>
iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::UniformBSplineCore ( const std::array< int64_t, parDim_ > &  ncoeffs,
const utils::TensorArray< geoDim_ > &  coeffs,
bool  clone = false,
Options< real_t >  options = Options<real_t>{} 
)
inline

Constructor for equidistant knot vectors.

Parameters
[in]ncoeffsNumber of coefficients per parametric dimension
[in]coeffsVectors of coefficients per parametric dimension
[in]cloneIf true, coefficients will be cloned. Otherwise, coefficients will be aliased
[in]optionsOptions configuration
Note
It is not checked whether vectors of coefficients are compatible with the given Options object if clone is false.

◆ UniformBSplineCore() [4/5]

template<typename real_t , short_t GeoDim, short_t... Degrees>
iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::UniformBSplineCore ( const std::array< int64_t, parDim_ > &  ncoeffs,
utils::TensorArray< geoDim_ > &&  coeffs,
Options< real_t >  options = Options<real_t>{} 
)
inline

Constructor for equidistant knot vectors.

Parameters
[in]ncoeffsNumber of coefficients per parametric dimension
[in]coeffsVectors of coefficients per parametric dimension
[in]optionsOptions configuration
Note
It is not checked whether vectors of coefficients are compatible with the given Options object if clone is false.

◆ UniformBSplineCore() [5/5]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename other_t >
iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::UniformBSplineCore ( const UniformBSplineCore< other_t, GeoDim, Degrees... > &  other,
Options< real_t >  options = Options<real_t>{} 
)
inline

Copy constructor.

Parameters
[in]otherUniform B-spline object to copy
[in]optionsOptions configuration

Member Function Documentation

◆ as_tensor()

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::Tensor iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::as_tensor ( ) const
inlineoverridevirtualnoexcept

Returns all coefficients as a single tensor.

Returns
Tensor of coefficients

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ as_tensor_()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<std::size_t... Is>
torch::Tensor iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::as_tensor_ ( std::index_sequence< Is... >  ) const
inlineprivatenoexcept

Returns all coefficients as a single tensor.

Returns
Tensor of coefficients

◆ as_tensor_size()

template<typename real_t , short_t GeoDim, short_t... Degrees>
int64_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::as_tensor_size ( ) const
inlineoverridevirtualnoexcept

Returns the size of the single tensor representation of all coefficients.

Returns
Size of the tensor

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ coeffs() [1/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
const utils::TensorArray< geoDim_ > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs ( ) const
inlinenoexcept

Returns a constant reference to the array of coefficient vectors.

Returns
Array of coefficient vectors

◆ coeffs() [2/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::TensorArray< geoDim_ > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs ( )
inlinenoexcept

Returns a non-constant reference to the array of coefficient vectors.

Returns
Array of coefficient vectord

◆ coeffs() [3/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
const torch::Tensor & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs ( short_t  i) const
inlinenoexcept

Returns a constant reference to the coefficient vector in the \(i\)-th dimension.

Parameters
[in]iGeometric dimension
Returns
Coefficient vector for the given geometric dimension

◆ coeffs() [4/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::Tensor & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs ( short_t  i)
inlinenoexcept

Returns a non-constant reference to the coefficient vector in the \(i\)-th dimension.

Parameters
[in]iGeometric dimension
Returns
Coefficient vector for the given geometric dimension

◆ coeffs_to_json()

template<typename real_t , short_t GeoDim, short_t... Degrees>
nlohmann::json iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs_to_json ( ) const
inline

Returns the B-spline object's coefficients as JSON object.

◆ coeffs_view() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::TensorArray< geoDim_ > iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs_view ( ) const
inlinenoexcept

Returns an array of views to the coefficient vectors.

Returns
Array of views to the coefficient vectors

◆ coeffs_view() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
const auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs_view ( short_t  i) const
inlinenoexcept

Returns a view to the coefficient vector in the \(i\)-th dimension.

Parameters
[in]iGeometric dimension
Returns
View of the coefficient vector for the given geometric dimension

◆ degree()

template<typename real_t , short_t GeoDim, short_t... Degrees>
static constexpr const short_t & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::degree ( short_t  i)
inlinestaticconstexprnoexcept

Returns a constant reference to the degree in the \(i\)-th dimension.

Parameters
[in]iParametric dimension
Returns
Degree for the given parametric dimension

◆ degrees()

template<typename real_t , short_t GeoDim, short_t... Degrees>
static constexpr const std::array< short_t, parDim_ > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::degrees ( )
inlinestaticconstexprnoexcept

Returns a constant reference to the array of degrees.

Returns
Array of degrees for all parametric dimensions

◆ device()

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::Device iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::device ( ) const
inlineoverridevirtualnoexcept

Returns the device property.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ device_index()

template<typename real_t , short_t GeoDim, short_t... Degrees>
int32_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::device_index ( ) const
inlineoverridevirtualnoexcept

Returns the device_index property.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ dtype()

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::Dtype iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::dtype ( ) const
inlineoverridevirtualnoexcept

Returns the dtype property.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ eval() [1/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval ( const torch::Tensor &  xi) const
inline

Returns the value of the B-spline object in the point xi

This implementation follows the procedure described in Chapters 2 and 3 of [2].

Algorithm: B-spline evaluation

  1. Determine the indices \((i_d)_{d=1}^{d_\text{par}}\) of the knot spans such that

\[ \boldsymbol{\xi} = \left(\xi_1, \dots, \xi_{d_\text{par}}\right)^\top \in \bigotimes_{d=1}^{d_\text{par}} [t_{i_d}, t_{i_d+1}). \]

  1. Evaluate the vectors of univariate B-spline basis functions (or their derivatives) that are non-zero at \(\boldsymbol{\xi}\)

\[ D^{r_d}\mathbf{B}_d = \left( D^{r_d} B_{i_d-p_d,p_d}, \dots, D^{r_d} B_{i_d,p_d} \right)^\top, \]

where \( p_d \) is the degree of the \(d\)-th univariate B-spline and \( r_d \) denotes the requested derivative in the \(d\)-direction.

  1. Multiply the tensor-product of the above row vectors by the column vector of control points

\[ \left( \bigotimes_{d=1}^{d_\text{par}} D^{r_d}\mathbf{B}_d \right) \cdot \mathbf{c}_\mathcal{J}, \]

where \(\mathcal{J}\) is the subset of global indices that belong to the coefficients

\[ \mathbf{c}_{i_1-p_1:i_1,\dots,i_\text{par}-p_\text{par}:i_\text{par}} \]

Template Parameters
derivComposition of derivative indicators of type deriv
Parameters
[in]xiPoint(s) where to evaluate the multivariate B-spline object
Returns
Value(s) of the multivariate B-spline evaluated at the point(s) xi

◆ eval() [2/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval ( const utils::TensorArray< parDim_ > &  xi) const
inline

Returns the value of the B-spline object in the point xi

This implementation follows the procedure described in Chapters 2 and 3 of [2].

Algorithm: B-spline evaluation

  1. Determine the indices \((i_d)_{d=1}^{d_\text{par}}\) of the knot spans such that

\[ \boldsymbol{\xi} = \left(\xi_1, \dots, \xi_{d_\text{par}}\right)^\top \in \bigotimes_{d=1}^{d_\text{par}} [t_{i_d}, t_{i_d+1}). \]

  1. Evaluate the vectors of univariate B-spline basis functions (or their derivatives) that are non-zero at \(\boldsymbol{\xi}\)

\[ D^{r_d}\mathbf{B}_d = \left( D^{r_d} B_{i_d-p_d,p_d}, \dots, D^{r_d} B_{i_d,p_d} \right)^\top, \]

where \( p_d \) is the degree of the \(d\)-th univariate B-spline and \( r_d \) denotes the requested derivative in the \(d\)-direction.

  1. Multiply the tensor-product of the above row vectors by the column vector of control points

\[ \left( \bigotimes_{d=1}^{d_\text{par}} D^{r_d}\mathbf{B}_d \right) \cdot \mathbf{c}_\mathcal{J}, \]

where \(\mathcal{J}\) is the subset of global indices that belong to the coefficients

\[ \mathbf{c}_{i_1-p_1:i_1,\dots,i_\text{par}-p_\text{par}:i_\text{par}} \]

Template Parameters
derivComposition of derivative indicators of type deriv
Parameters
[in]xiPoint(s) where to evaluate the multivariate B-spline object
Returns
Value(s) of the multivariate B-spline evaluated at the point(s) xi

◆ eval() [3/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval ( const utils::TensorArray< parDim_ > &  xi,
const utils::TensorArray< parDim_ > &  knot_indices 
) const
inline

Returns the value of the univariate B-spline object in the points xi

This function implements steps 2-3 of algorithm BSplineEvaluation for univariate B-splines (i.e. \(d_\text{par}=1\))

Template Parameters
derivComposition of derivative indicators of type deriv
Parameters
[in]xiPoint(s) where to evaluate the univariate B-spline object
[in]knot_indicesKnot indices where to evaluate the univariate B-spline object
Returns
Value(s) of the univariate B-spline evaluated at the point(s) xi

◆ eval() [4/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval ( const utils::TensorArray< parDim_ > &  xi,
const utils::TensorArray< parDim_ > &  knot_indices,
const torch::Tensor &  coeff_indices 
) const
inline

Returns the value of the univariate B-spline object in the points xi

This function implements steps 2-3 of algorithm BSplineEvaluation for univariate B-splines (i.e. \(d_\text{par}=1\))

Template Parameters
derivComposition of derivative indicators of type deriv
Parameters
[in]xiPoint(s) where to evaluate the univariate B-spline object
[in]knot_indicesKnot indices where to evaluate the univariate B-spline object
[in]coeff_indicesCoefficient indices where to evaluate the univariate B-spline object
Returns
Value(s) of the univariate B-spline evaluated at the point(s) xi

◆ eval_basfunc() [1/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_basfunc ( const torch::Tensor &  xi) const
inline

Returns the vector of multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

Parameters
[in]xiPoint(s) where to evaluate the B-spline object
Returns
Multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

◆ eval_basfunc() [2/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_basfunc ( const torch::Tensor &  xi,
const torch::Tensor &  knot_indices 
) const
inline

Returns the vector of multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

Parameters
[in]xiPoint(s) where to evaluate the B-spline object
[in]knot_indicesKnot indices where to evaluate the univariate B-spline object
Returns
Multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

◆ eval_basfunc() [3/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_basfunc ( const utils::TensorArray< parDim_ > &  xi) const
inline

Returns the vector of multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

Parameters
[in]xiPoint(s) where to evaluate the B-spline object
Returns
Multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

◆ eval_basfunc() [4/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<deriv deriv = deriv::func, bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_basfunc ( const utils::TensorArray< parDim_ > &  xi,
const utils::TensorArray< parDim_ > &  knot_indices 
) const
inline

Returns the vector of multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

Parameters
[in]xiPoint(s) where to evaluate the B-spline object
[in]knot_indicesKnot indices where to evaluate the univariate B-spline object
Returns
Multivariate B-spline basis functions (or their derivatives) evaluated in the point xi

◆ eval_basfunc_univariate()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<short_t degree, short_t dim, short_t deriv>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_basfunc_univariate ( const torch::Tensor &  xi,
const torch::Tensor &  knot_indices 
) const
inlineprotected

Returns the vector of univariate B-spline basis functions (or their derivatives) evaluated in the point xi

This function implements step 2 of algorithm BSplineEvaluation, that is, it evaluates the vector of univariate B-spline basis functions (or their derivatives) that are non-zero at \(\xi_d \in [t_{i_d}, t_{i_d+1})\)

\[ D^{r_d}\mathbf{B}_d(\xi_d) = \left( D^{r_d} B_{i_d-p_d,p_d}(\xi_d), \dots, D^{r_d} B_{i_d,p_d}(\xi_d) \right)^\top, \]

where \( p_d \) is the degree of the \(d\)-th univariate B-spline and \( r_d \) denotes the requested derivative in the \(d\)-direction.

According to the procedure described in Chapters 2 and 3 of [2] this can be accomplished by the following expression

\[ D^{r_d}\mathbf{B}_d(\xi_d) = \frac{p_d!}{(p_d-r_d)!}\mathbf{R}_1(\xi_d)\cdot \cdots \cdot \mathbf{R}_{p_d-r_d}(\xi_d) D\mathbf{R}_{p_d-r_d+1}\cdot \cdots \cdot D\mathbf{R}_{p_d}(\xi_d), \]

where (cf. Equation (2.20) in [2])

\[ \mathbf{R}_k(\xi_d) = \begin{pmatrix} \frac{t_{i_p+1} - \xi_d}{t_{i_p+1} - t_{i_p+1-k}} & \frac{\xi_d - t_{i_p+1-k}}{t_{i_p+1} - t_{i_p+1-k}} & 0 & \cdots & 0 \\ 0 & \frac{t_{i_p+2} - \xi_d}{t_{i_p+2} - t_{i_p+2-k}} & \frac{\xi_d - t_{i_p+2-k}}{t_{i_p+2} - t_{i_p+1-k}} & \cdots & 0 \\ \vdots & \vdots & \ddots & \ddots & \vdots \\ 0 & 0 & \cdots & \frac{t_{i_p+k} - \xi_d}{t_{i_p+k} - t_{i_p}} & \frac{\xi_d - t_{i_p}}{t_{i_p+k} - t_{i_p}} \end{pmatrix} \]

and (cf. Equation (3.30) in [2])

\[ D\mathbf{R}_k(\xi_d) = \begin{pmatrix} \frac{-1}{t_{i_p+1} - t_{i_p+1-k}} & \frac{1}{t_{i_p+1} - t_{i_p+1-k}} & 0 & \cdots & 0 \\ 0 & \frac{-1}{t_{i_p+2} - t_{i_p+2-k}} & \frac{1}{t_{i_p+2} - t_{i_p+1-k}} & \cdots & 0 \\ \vdots & \vdots & \ddots & \ddots & \vdots \\ 0 & 0 & \cdots & \frac{-1}{t_{i_p+k} - t_{i_p}} & \frac{1}{t_{i_p+k} - t_{i_p}} \end{pmatrix}. \]

To improve computational efficiency, the prefactor

\[ \frac{p_d!}{(p_d-r_d)!}=p_d \cdots (p_d-r_d+1) \]

is computed as compile-time expression by the eval_prefactor() function.

Moreover, the above expression for \(D^{r_d}\mathbf{B}_d(\xi_d)\) is evaluated as described in Algorithm 2.22 (R-vector version) in [2]) and its generalization to derivatives, respectively.

The algorithm goes as follows:

  1. \(\mathbf{b} = 1\)
  2. For \(k = 1, \dots, p_d-r_d\)
    1. \(\mathbf{t}_1 = \left(t_{i_d-k+1},\dots,t_{i_d}\right)\)
    2. \(\mathbf{t}_2 = \left(t_{i_d+1},\dots,t_{i_d+k}\right)\)
    3. \(\mathbf{w} = \left(\xi_d-\mathbf{t}_1\right)\div\left(\mathbf{t}_2-\mathbf{t}_1\right)\)
    4. \(\mathbf{b} = \left[\left(1-\mathbf{w}\right)\odot\mathbf{b}, 0\right] + \left[0, \mathbf{w}\odot\mathbf{b}\right]\)
  3. For \(k = p_d-r_d+1, \dots, p_d\)
    1. \(\mathbf{t}_1 = \left(t_{i_d-k+1},\dots,t_{i_d}\right)\)
    2. \(\mathbf{t}_2 = \left(t_{i_d+1},\dots,t_{i_d+k}\right)\)
    3. \(\mathbf{w} = 1\div\left(\mathbf{t}_2-\mathbf{t}_1\right)\)
    4. \(\mathbf{b} = \left[-\mathbf{w}\odot\mathbf{b}, 0\right] + \left[0, \mathbf{w}\odot\mathbf{b}\right]\)

where \(\div\) and \(\odot\) denote the element-wise division and multiplication of vectors, respectively.

◆ eval_from_precomputed() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::BlockTensor< torch::Tensor, 1, geoDim_ > iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_from_precomputed ( const torch::Tensor &  basfunc,
const torch::Tensor &  coeff_indices,
int64_t  numeval,
torch::IntArrayRef  sizes 
) const
inlineoverridevirtual

Returns the value of the B-spline object from precomputed basis function.

This function implements steps 2-3 of algorithm BSplineEvaluation for univariate B-splines (i.e. \(d_\text{par}=1\))

Parameters
[in]basfuncValue(s) of the multivariate B-spline basis functions evaluated at the point(s) xi
[in]coeff_indicesIndices where to evaluate the coefficients
[in]numevalNumber of evaluation points
[in]sizesDimension of the result
Returns
Value(s) of the univariate B-spline object
Note
This function does not work of the basis functions are evaluated with memory_optimized flag to true

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ eval_from_precomputed() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::BlockTensor< torch::Tensor, 1, geoDim_ > iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_from_precomputed ( const utils::TensorArray< parDim_ > &  basfunc,
const torch::Tensor &  coeff_indices,
int64_t  numeval,
torch::IntArrayRef  sizes 
) const
inlineoverride

Returns the value of the B-spline object from precomputed basis function.

This function implements steps 2-3 of algorithm BSplineEvaluation for univariate B-splines (i.e. \(d_\text{par}=1\))

Parameters
[in]basfuncValue(s) of the multivariate B-spline basis functions evaluated at the point(s) xi
[in]coeff_indicesIndices where to evaluate the coefficients
[in]numevalNumber of evaluation points
[in]sizesDimension of the result
Returns
Value(s) of the univariate B-spline object
Note
This function does not work of the basis functions are evaluated with memory_optimized flag to true

◆ eval_prefactor()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<int64_t degree, int64_t deriv, int64_t terminal = degree - deriv>
int64_t constexpr iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::eval_prefactor ( ) const
inlineconstexprprivate

Computes the prefactor \(p_d!/(p_d-r_d)! = p_d \cdots (p_d-r_d+1)\).

◆ find_coeff_indices() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::find_coeff_indices ( const torch::Tensor &  indices) const
inline

Returns the indices of the coefficients corresponding to the knot indices indices

Parameters
[in]indicesIndices of the knot spans
Returns
Indices of the coefficients corresponding to the knot indices

◆ find_coeff_indices() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<bool memory_optimized = false>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::find_coeff_indices ( const utils::TensorArray< parDim_ > &  indices) const
inline

Returns the indices of the coefficients corresponding to the knot indices indices

Parameters
[in]indicesIndices of the knot spans
Returns
Indices of the coefficients corresponding to the knot indices

◆ find_knot_indices() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::find_knot_indices ( const torch::Tensor &  xi) const
inlinenoexcept

Returns the indices of knot spans containing xi

This function returns the indices \((i_d)_{d=1}^{d_\text{par}}\) of the knot spans such that

\[ \boldsymbol{\xi} \in [t_{i_1}, t_{i_1+1}) \times [t_{i_2}, t_{i_2+1}) \times \dots \times [t_{i_{d_\text{par}}}, t_{i_{d_\text{par}}+1}). \]

The indices are returned as utils::TensorArray<parDim_> in the same order as provided in xi

Parameters
[in]xiPoint(s) where to evaluate the B-spline object
Returns
Indices of the knot spans containing xi

◆ find_knot_indices() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::TensorArray< parDim_ > iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::find_knot_indices ( const utils::TensorArray< parDim_ > &  xi) const
inlinenoexcept

Returns the indices of knot spans containing xi

This function returns the indices \((i_d)_{d=1}^{d_\text{par}}\) of the knot spans such that

\[ \boldsymbol{\xi} \in [t_{i_1}, t_{i_1+1}) \times [t_{i_2}, t_{i_2+1}) \times \dots \times [t_{i_{d_\text{par}}}, t_{i_{d_\text{par}}+1}). \]

The indices are returned as utils::TensorArray<parDim_> in the same order as provided in xi

Parameters
[in]xiPoint(s) where to evaluate the B-spline object
Returns
Indices of the knot spans containing xi

◆ from_gismo()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename BSpline >
auto & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::from_gismo ( BSpline bspline,
bool  ,
bool   
)
inline

◆ from_json()

template<typename real_t , short_t GeoDim, short_t... Degrees>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::from_json ( const nlohmann::json &  json)
inline

Updates the B-spline object from JSON object.

◆ from_tensor()

template<typename real_t , short_t GeoDim, short_t... Degrees>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::from_tensor ( const torch::Tensor &  tensor)
inlineoverridevirtualnoexcept

Sets all coefficients from a single tensor.

Parameters
[in]tensorTensor from which to extract the coefficients
Returns
Updated spline object

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ from_tensor_()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<std::size_t... Is>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::from_tensor_ ( std::index_sequence< Is... >  ,
const torch::Tensor &  tensor 
)
inlineprivatenoexcept

Sets all coefficients from a single tensor.

Returns
Updates spline object

◆ from_xml() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::from_xml ( const pugi::xml_document &  doc,
int  id = 0,
std::string  label = "",
int  index = -1 
)
inline

Updates the B-spline object from XML object.

◆ from_xml() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::from_xml ( const pugi::xml_node &  root,
int  id = 0,
std::string  label = "",
int  index = -1 
)
inline

Updates the B-spline object from XML node.

◆ geoDim()

template<typename real_t , short_t GeoDim, short_t... Degrees>
static constexpr short_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::geoDim ( )
inlinestaticconstexprnoexcept

Returns the geometric dimension.

Returns
Number of geometric dimensions

◆ greville()

template<typename real_t , short_t GeoDim, short_t... Degrees>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::greville ( bool  interior = false) const
inline

Returns the Greville abscissae.

The Greville abscissae are defined as

\[ g_{i_d} = \frac{\xi_{i_d+1} + \xi_{i_d+2} + \dots + \xi_{i_d+p_d+1}}{p_d-1} \]

Parameters
[in]interiorIf true only interior Greville abscissae are considered
Returns
Array of Greville abscissae

◆ init_coeffs()

template<typename real_t , short_t GeoDim, short_t... Degrees>
void iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::init_coeffs ( enum init  init)
inline

Initializes the B-spline coefficients.

◆ init_knots()

template<typename real_t , short_t GeoDim, short_t... Degrees>
void iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::init_knots ( )
inline

Initializes the B-spline knots.

◆ is_nonuniform()

template<typename real_t , short_t GeoDim, short_t... Degrees>
static constexpr bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::is_nonuniform ( )
inlinestaticconstexprnoexcept

Returns true if the B-spline is non-uniform.

◆ is_sparse()

template<typename real_t , short_t GeoDim, short_t... Degrees>
bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::is_sparse ( ) const
inlineoverridevirtualnoexcept

Returns true if the layout is sparse.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ is_uniform()

template<typename real_t , short_t GeoDim, short_t... Degrees>
static constexpr bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::is_uniform ( )
inlinestaticconstexprnoexcept

Returns true if the B-spline is uniform.

◆ isclose()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename other_t , short_t GeoDim_, short_t... Degrees_>
bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::isclose ( const UniformBSplineCore< other_t, GeoDim_, Degrees_... > &  other,
real_t  rtol = real_t{1e-5},
real_t  atol = real_t{1e-8} 
) const
inline

Returns true if both B-spline objects are close up to the given tolerances.

◆ knots() [1/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
const utils::TensorArray< parDim_ > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::knots ( ) const
inlinenoexcept

Returns a constant reference to the array of knot vectors.

Returns
Array of knot vectors

◆ knots() [2/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::TensorArray< parDim_ > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::knots ( )
inlinenoexcept

Returns a non-constant reference to the array of knot vectors.

Returns
Array of knot vectors

◆ knots() [3/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
const torch::Tensor & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::knots ( short_t  i) const
inlinenoexcept

Returns a constant reference to the knot vector in the \(i\)-th dimension.

Parameters
[in]iParametric dimension
Returns
Knot vector for the given parametric dimension

◆ knots() [4/4]

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::Tensor & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::knots ( short_t  i)
inlinenoexcept

Returns a non-constant reference to the knot vector in the \(i\)-th dimension.

Parameters
[in]iParametric dimension
Returns
Knot vector for the given parametric dimension

◆ knots_to_json()

template<typename real_t , short_t GeoDim, short_t... Degrees>
nlohmann::json iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::knots_to_json ( ) const
inline

Returns the B-spline object's knots as JSON object.

◆ layout()

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::Layout iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::layout ( ) const
inlineoverridevirtualnoexcept

Returns the layout property.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ load()

template<typename real_t , short_t GeoDim, short_t... Degrees>
void iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::load ( const std::string &  filename,
const std::string &  key = "bspline" 
)
inline

Loads the B-spline from file.

◆ ncoeffs() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
const std::array< int64_t, parDim_ > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::ncoeffs ( ) const
inlinenoexcept

Returns a constant reference to the array of coefficient vector dimensions.

Returns
Array of coefficient vector dimensions

◆ ncoeffs() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
int64_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::ncoeffs ( short_t  i) const
inlinenoexcept

Returns the total number of coefficients in the \(i\)-th direction.

Parameters
[in]iParametric dimension
Returns
Total number of coefficients in given parametric dimension

◆ ncumcoeffs()

template<typename real_t , short_t GeoDim, short_t... Degrees>
int64_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::ncumcoeffs ( ) const
inlinenoexcept

Returns the total number of coefficients.

Returns
Total number of coefficients

◆ nknots() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
const std::array< int64_t, parDim_ > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::nknots ( ) const
inlinenoexcept

Returns a constant reference to the array of knot vector dimensions.

Returns
Array of knot vector dimensions

◆ nknots() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
int64_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::nknots ( short_t  i) const
inlinenoexcept

Returns the dimension of the knot vector in the \(i\)-th dimension.

Parameters
[in]iParametric dimension
Returns
Knot vector dimension for the given parametric dimension

◆ operator!=()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename other_t , short_t GeoDim_, short_t... Degrees_>
bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::operator!= ( const UniformBSplineCore< other_t, GeoDim_, Degrees_... > &  other) const
inline

Returns true if both B-spline objects are different.

◆ operator==()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename other_t , short_t GeoDim_, short_t... Degrees_>
bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::operator== ( const UniformBSplineCore< other_t, GeoDim_, Degrees_... > &  other) const
inline

Returns true if both B-spline objects are the same.

◆ options()

template<typename real_t , short_t GeoDim, short_t... Degrees>
const Options< real_t > & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::options ( ) const
inlinenoexcept

Returns a constant reference to the B-spline object's options.

◆ parDim()

template<typename real_t , short_t GeoDim, short_t... Degrees>
static constexpr short_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::parDim ( )
inlinestaticconstexprnoexcept

Returns the parametric dimension.

Returns
Number of parametric dimensions

◆ pinned_memory()

template<typename real_t , short_t GeoDim, short_t... Degrees>
bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::pinned_memory ( ) const
inlineoverridevirtualnoexcept

Returns the pinned_memory property.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ read()

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::serialize::InputArchive & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::read ( torch::serialize::InputArchive &  archive,
const std::string &  key = "bspline" 
)
inline

Reads the B-spline from a torch::serialize::InputArchive object.

◆ requires_grad()

template<typename real_t , short_t GeoDim, short_t... Degrees>
bool iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::requires_grad ( ) const
inlineoverridevirtualnoexcept

Returns the requires_grad property.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ save()

template<typename real_t , short_t GeoDim, short_t... Degrees>
void iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::save ( const std::string &  filename,
const std::string &  key = "bspline" 
) const
inline

Saves the B-spline to file.

◆ set_requires_grad()

template<typename real_t , short_t GeoDim, short_t... Degrees>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::set_requires_grad ( bool  requires_grad)
inlineoverridevirtualnoexcept

Sets the B-spline object's requires_grad property.

Note
: It is only necessary to set requires_grad to true if gradients with respect to B-spline entities, e.g., the control points should be computed. For computing the gradients with respect to the sampling points the B-spline's requires_grad property can be false.

Implements iganet::BSplinePatch< real_t, GeoDim, sizeof...(Degrees)>.

◆ to_gismo() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::to_gismo ( ) const
inline

Converts the B-spline object into a gsBSpline object of the parametric dimension is one and a gsTensorBSpline object otherwise.

◆ to_gismo() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename BSpline >
BSpline & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::to_gismo ( BSpline bspline,
bool  ,
bool   
) const
inline

◆ to_json()

template<typename real_t , short_t GeoDim, short_t... Degrees>
nlohmann::json iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::to_json ( ) const
inlineoverridevirtual

Returns the B-spline object as JSON object.

Implements iganet::utils::Serializable.

◆ to_xml() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
pugi::xml_document iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::to_xml ( int  id = 0,
std::string  label = "",
int  index = -1 
) const
inline

Returns the B-spline object as XML object.

◆ to_xml() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
pugi::xml_node & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::to_xml ( pugi::xml_node &  root,
int  id = 0,
std::string  label = "",
int  index = -1 
) const
inline

Returns the B-spline object as XML node.

◆ transform() [1/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::transform ( const std::function< std::array< real_t, geoDim_ >(const std::array< real_t, parDim_ > &)>  transformation)
inline

Transforms the coefficients based on the given mapping.

◆ transform() [2/2]

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<std::size_t N>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::transform ( const std::function< std::array< real_t, N >(const std::array< real_t, parDim_ > &)>  transformation,
std::array< short_t, N dims 
)
inline

Transforms the coefficients based on the given mapping.

◆ uniform_refine()

template<typename real_t , short_t GeoDim, short_t... Degrees>
UniformBSplineCore & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::uniform_refine ( int  numRefine = 1,
int  dim = -1 
)
inline

Returns the B-spline object with uniformly refined knot and coefficient vectors.

If dim = -1, new knot values are inserted uniformly in each knot span in all spatial dimensions. Otherwise, i.e., dim != -1 new knots are only inserted in the specified dimension.

◆ update_coeffs()

template<typename real_t , short_t GeoDim, short_t... Degrees>
void iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::update_coeffs ( const utils::TensorArray< parDim_ > &  knots,
const utils::TensorArray< parDim_ > &  knot_indices 
)
inlineprotected

Updates the B-spline coefficients after knot insertion.

◆ update_coeffs_univariate()

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<short_t degree, short_t dim>
auto iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::update_coeffs_univariate ( const torch::Tensor &  knots,
const torch::Tensor &  knot_indices 
) const
inlineprotected

Returns the knot insertion matrix.

This functions implements the Oslo algorithm (Algorithm 4.11 in [2]) to compute the univariate knot insertion matrix from the given knot vector to the new knot vector passed as argument knots.

◆ write()

template<typename real_t , short_t GeoDim, short_t... Degrees>
torch::serialize::OutputArchive & iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::write ( torch::serialize::OutputArchive &  archive,
const std::string &  key = "bspline" 
) const
inline

Writes the B-spline into a torch::serialize::OutputArchive object.

Friends And Related Symbol Documentation

◆ BSplineCommon

template<typename real_t , short_t GeoDim, short_t... Degrees>
template<typename BSplineCore >
friend class BSplineCommon
friend

Enable access to private members.

Member Data Documentation

◆ coeffs_

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::TensorArray<geoDim_> iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::coeffs_
protected

Array storing the coefficients of the control net \(\left(\mathbf{c}_{i_d}\right)_{i_d=1}^{n_d}\), \(\mathbf{c}_{i_d}\in\mathbb{R}^{d_\text{geo}}\).

◆ degrees_

template<typename real_t , short_t GeoDim, short_t... Degrees>
constexpr const std::array<short_t, parDim_> iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::degrees_ = {Degrees...}
staticconstexprprotected

Array storing the degrees \(\left(p_d\right)_{d=1}^{d_\text{par}}\).

◆ geoDim_

template<typename real_t , short_t GeoDim, short_t... Degrees>
constexpr const short_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::geoDim_ = GeoDim
staticconstexprprotected

Dimension of the geometric space \(\Omega\subset\mathbb{R}^{d_\text{geo}}\).

◆ knots_

template<typename real_t , short_t GeoDim, short_t... Degrees>
utils::TensorArray<parDim_> iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::knots_
protected

Array storing the knot vectors \(\left(\left(t_{i_d}\right)_{i_d=1}^{n_d+p_d+1}\right)_{d=1}^{d_\text{par}}\).

◆ ncoeffs_

template<typename real_t , short_t GeoDim, short_t... Degrees>
std::array<int64_t, parDim_> iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::ncoeffs_
protected

Array storing the sizes of the coefficients of the control net \(\left(n_d\right)_{d=1}^{d_\text{par}}\).

◆ ncoeffs_reverse_

template<typename real_t , short_t GeoDim, short_t... Degrees>
std::array<int64_t, parDim_> iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::ncoeffs_reverse_
protected

Array storing the sizes of the coefficients of the control net \(\left(n_d\right)_{d=1}^{d_\text{par}}\) in reverse order (needed for coeffs_view)

◆ nknots_

template<typename real_t , short_t GeoDim, short_t... Degrees>
std::array<int64_t, parDim_> iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::nknots_
protected

Array storing the sizes of the knot vectors \(\left(n_d+p_d+1\right)_{d=1}^{d_\text{par}}\).

◆ options_

template<typename real_t , short_t GeoDim, short_t... Degrees>
Options<real_t> iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::options_
protected

◆ parDim_

template<typename real_t , short_t GeoDim, short_t... Degrees>
constexpr const short_t iganet::UniformBSplineCore< real_t, GeoDim, Degrees >::parDim_ = sizeof...(Degrees)
staticconstexprprotected

Dimension of the parametric space \(\hat\Omega=[0,1]^{d_\text{par}}\).


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