IgANet
IgANets - Isogeometric Analysis Networks
Loading...
Searching...
No Matches
nurbs.hpp
Go to the documentation of this file.
1
15#pragma once
16
17#include <bspline.hpp>
18
19namespace iganet {
20
33//
39//
48template <typename real_t, short_t GeoDim, short_t... Degrees>
50 : public UniformBSplineCore<real_t, GeoDim + 1, Degrees...> {
51private:
53 using Base = UniformBSplineCore<real_t, GeoDim + 1, Degrees...>;
54
55public:
57 using value_type = real_t;
58
64 template <template <typename, short_t, short_t...> class BSpline,
65 std::make_signed<short_t>::type degree_elevate = 0>
67 BSpline<real_t, GeoDim + 1, (Degrees + degree_elevate)...>;
68
72 using self_type =
74
78 template <typename other_t, short_t GeoDim_, short_t... Degrees_>
80
83 template <typename other_t>
85
87 static constexpr bool is_uniform() { return true; }
88
90 static constexpr bool is_nonuniform() { return false; }
91
97
105 UniformNurbsCore(const std::array<int64_t, Base::parDim_> &ncoeffs,
106 enum init init = init::greville,
108 : Base(ncoeffs, init, options) {
110 Base::coeffs_[GeoDim] = torch::ones_like(Base::coeffs_[GeoDim]);
111 }
112
126 UniformNurbsCore(const std::array<int64_t, Base::parDim_> &ncoeffs,
128 bool clone = false,
130 : Base(ncoeffs, coeffs, clone, options) {}
131
146
152 template <typename other_t>
158
164 inline static constexpr short_t geoDim() noexcept { return GeoDim; }
165
171 inline const torch::Tensor &weights() const noexcept {
172 return Base::coeffs_[GeoDim];
173 }
174
180 inline torch::Tensor &weights() noexcept { return Base::coeffs_[GeoDim]; }
181};
182
189template <typename real_t, short_t GeoDim, short_t... Degrees>
191 : public NonUniformBSplineCore<real_t, GeoDim, Degrees...> {
192private:
195
196protected:
198 torch::Tensor weights_;
199
200public:
202 using value_type = real_t;
203
209 template <template <typename, short_t, short_t...> class BSpline,
210 std::make_signed<short_t>::type degree_elevate = 0>
212
216 using self_type =
218
222 template <typename other_t, short_t GeoDim_, short_t... Degrees_>
224
227 template <typename other_t>
230
232 static constexpr bool is_uniform() { return false; }
233
235 static constexpr bool is_nonuniform() { return true; }
236};
237
239template <typename real_t, short_t GeoDim, short_t... Degrees>
242
244template <typename real_t, short_t GeoDim, short_t... Degrees>
245inline std::ostream &
246operator<<(std::ostream &os,
249 return os;
250}
251
253template <typename real_t, short_t GeoDim, short_t... Degrees>
256
258template <typename real_t, short_t GeoDim, short_t... Degrees>
259inline std::ostream &
260operator<<(std::ostream &os,
263 return os;
264}
265} // namespace iganet
Multivariate B-splines.
B-spline (common high-level functionality)
Definition bspline.hpp:3261
Tensor-product non-uniform B-spline (core functionality)
Definition bspline.hpp:2720
Tensor-product non-uniform rational B-spline with non-uniform knot vectors (core functionality)
Definition nurbs.hpp:191
static constexpr bool is_nonuniform()
Returns true if the B-spline is non-uniform.
Definition nurbs.hpp:235
torch::Tensor weights_
Tensor storing the weights.
Definition nurbs.hpp:198
BSpline< real_t, GeoDim,(Degrees+degree_elevate)... > derived_type
Deduces the type of the template template parameter BSpline when exposed to the class template parame...
Definition nurbs.hpp:211
real_t value_type
Value type.
Definition nurbs.hpp:202
static constexpr bool is_uniform()
Returns true if the B-spline is uniform.
Definition nurbs.hpp:232
typename Base::template derived_type< NonUniformNurbsCore, degree_elevate > self_type
Deduces the self-type possibly degrees (de-)elevated by the additive constant degree_elevate
Definition nurbs.hpp:217
The Options class handles the automated determination of dtype from the template argument and the sel...
Definition options.hpp:90
Tensor-product uniform B-spline (core functionality)
Definition bspline.hpp:182
const Options< real_t > & options() const noexcept
Returns a constant reference to the B-spline object's options.
Definition bspline.hpp:321
utils::TensorArray< geoDim_ > coeffs_
Array storing the coefficients of the control net , .
Definition bspline.hpp:219
const utils::TensorArray< geoDim_ > & coeffs() const noexcept
Returns a constant reference to the array of coefficient vectors.
Definition bspline.hpp:534
const std::array< int64_t, parDim_ > & ncoeffs() const noexcept
Returns a constant reference to the array of coefficient vector dimensions.
Definition bspline.hpp:609
Tensor-product non-uniform rational B-spline with uniform knot vector (core functionality)
Definition nurbs.hpp:50
const torch::Tensor & weights() const noexcept
Returns a constant reference to the weights.
Definition nurbs.hpp:171
UniformBSplineCore< real_t, GeoDim+1, Degrees... > Base
Base type.
Definition nurbs.hpp:53
typename Base::template derived_type< UniformNurbsCore, degree_elevate > self_type
Deduces the self-type possibly degrees (de-)elevated by the additive constant degree_elevate
Definition nurbs.hpp:73
static constexpr short_t geoDim() noexcept
Definition nurbs.hpp:164
UniformNurbsCore(Options< real_t > options=Options< real_t >{})
Default constructor.
Definition nurbs.hpp:95
UniformNurbsCore(const std::array< int64_t, Base::parDim_ > &ncoeffs, const utils::TensorArray< Base::geoDim_ > &coeffs, bool clone=false, Options< real_t > options=Options< real_t >{})
Constructor for equidistant knot vectors.
Definition nurbs.hpp:126
static constexpr bool is_uniform()
Returns true if the B-spline is uniform.
Definition nurbs.hpp:87
UniformNurbsCore(const UniformNurbsCore< other_t, GeoDim+1, Degrees... > &other, Options< real_t > options=Options< real_t >{})
Copy constructor.
Definition nurbs.hpp:153
UniformNurbsCore(const std::array< int64_t, Base::parDim_ > &ncoeffs, utils::TensorArray< Base::geoDim_ > &&coeffs, Options< real_t > options=Options< real_t >{})
Constructor for equidistant knot vectors.
Definition nurbs.hpp:142
BSpline< real_t, GeoDim+1,(Degrees+degree_elevate)... > derived_type
Deduces the type of the template template parameter BSpline when exposed to the class template parame...
Definition nurbs.hpp:67
static constexpr bool is_nonuniform()
Returns true if the B-spline is non-uniform.
Definition nurbs.hpp:90
UniformNurbsCore(const std::array< int64_t, Base::parDim_ > &ncoeffs, enum init init=init::greville, Options< real_t > options=Options< real_t >{})
Constructor for equidistant knot vectors.
Definition nurbs.hpp:105
real_t value_type
Value type.
Definition nurbs.hpp:57
torch::Tensor & weights() noexcept
Returns a non-constant reference to the weights.
Definition nurbs.hpp:180
std::array< torch::Tensor, N > TensorArray
Definition tensorarray.hpp:28
Definition boundary.hpp:22
constexpr bool is_SplineType_v
Alias to the value of is_SplineType.
Definition bspline.hpp:3243
std::ostream & operator<<(std::ostream &os, const Boundary< Spline > &obj)
Print (as string) a Boundary object.
Definition boundary.hpp:1978
init
Enumerator for specifying the initialization of B-spline coefficients.
Definition bspline.hpp:55
short int short_t
Definition core.hpp:74
virtual void pretty_print(std::ostream &os=Log(log::info)) const =0
Returns a string representation of the object.