Medusa  1.1
Coordinate Free Mehless Method implementation
NURBSPatch.hpp
Go to the documentation of this file.
1 #ifndef MEDUSA_BITS_DOMAINS_NURBSPATCH_HPP_
2 #define MEDUSA_BITS_DOMAINS_NURBSPATCH_HPP_
3 
9 #include "NURBSPatch_fwd.hpp"
10 
13 
14 namespace mm {
15 
16 template <typename vec_t, typename param_vec_t>
18  const std::array<Range<scalar_t>, param_dim>& ks,
19  const std::array<int, param_dim>& in_p) :
20  weighted_control_points{ wcp }, knots{ ks }, p{ in_p },
21  der_structure{ nullptr } {
22  assert_msg(param_dim <= 2, "Only NURBS curves and surfaces are currently supported.");
23 }
24 
25 template <typename vec_t, typename param_vec_t>
27  const std::array<Range<scalar_t>, param_dim>& ks,
28  const std::array<int, param_dim>& in_p) :
29  knots{ ks }, p{ in_p }, der_structure{ nullptr } {
30  assert_msg(param_dim <= 2, "Only NURBS curves and surfaces are currently supported.");
31 
32  weighted_control_points = nurbs_patch_internal::NURBSPatchHelper<vec_t, param_vec_t>::
33  cpToWcp(cp, w);
34 }
35 
36 template <typename vec_t, typename param_vec_t>
38  : weighted_control_points{ patch.weighted_control_points },
39  knots{ patch.knots }, p { patch.p } {
40  if (patch.der_structure != nullptr) {
41  der_structure = std::unique_ptr<std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>>(
42  new std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>(*patch.der_structure));
43  }
44 }
45 
46 template <typename vec_t, typename param_vec_t>
48  const NURBSPatch<vec_t, param_vec_t>& other) noexcept {
49  weighted_control_points = other.weighted_control_points;
50  p = other.p;
51  knots = other.knots;
52 
53  if (other.der_structure == nullptr && der_structure != nullptr) {
54  der_structure = nullptr;
55  } else if (other.der_structure != nullptr) {
56  der_structure = std::unique_ptr<std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>>(
57  new std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>(*other.der_structure));
58  }
59 
60  return *this;
61 }
62 
63 template <typename vec_t, typename param_vec_t>
66 }
67 
68 template <typename vec_t, typename param_vec_t>
69 vec_t NURBSPatch<vec_t, param_vec_t>::evaluate(const param_vec_t& t, scalar_t* w) const {
70  proj_vec_t weighted_pt = evaluateWeighted(t);
71  // Project into the desired dimension.
72  vec_t pt;
73  for (int i = 0; i < dim; i++) {
74  pt(i) = weighted_pt(i) / weighted_pt(proj_dim - 1);
75  }
76  *w = weighted_pt(proj_dim - 1);
77  return pt;
78 }
79 
80 template <typename vec_t, typename param_vec_t>
81 vec_t NURBSPatch<vec_t, param_vec_t>::evaluate(const param_vec_t& t) const {
82  scalar_t w;
83  return evaluate(t, &w);
84 }
85 
86 template <typename vec_t, typename param_vec_t>
90 }
91 
92 template <typename vec_t, typename param_vec_t>
95 NURBSPatch<vec_t, param_vec_t>::jacobian(const param_vec_t& t, const vec_t& pt,
96  const scalar_t& w) const {
98 }
99 
100 template <typename vec_t, typename param_vec_t>
103 NURBSPatch<vec_t, param_vec_t>::jacobian(const param_vec_t& t) const {
104  scalar_t w;
105  vec_t pt = evaluate(t, &w);
106  return jacobian(t, pt, w);
107 }
108 
109 template <typename vec_t, typename param_vec_t>
110 std::pair<vec_t,
115  scalar_t w;
116  vec_t pt = evaluate(t, &w);
117  Eigen::Matrix<scalar_t, dim, param_dim> jm = jacobian(t, pt, w);
118  return std::pair<vec_t, Eigen::Matrix<scalar_t, dim, param_dim>>(pt, jm);
119 }
120 
121 template <typename vec_t, typename param_vec_t>
123  param_vec_t beg, end;
124  for (int i = 0; i < param_dim; i++) {
125  beg(i) = knots[i].front();
126  end(i) = knots[i].back();
127  }
128 
129  return BoxShape<param_vec_t>(beg, end);
130 }
131 
132 template <typename vec_t, typename param_vec_t>
137 }
138 
139 template <typename vec_t, typename param_vec_t>
141  const Vec<scalar_t, param_dim - 1>& t, int i, scalar_t epsilon) const {
143  getPatchParameterFromBoundaryParameter(*this, t, i, epsilon);
144 }
145 
146 template <typename vec_t, typename param_vec_t>
149  getBoundaryParameterFromPatchParameter(const param_vec_t& t, int i) const {
152 }
153 
155 template <typename vec_t>
156 struct nurbs_patch_internal::NURBSPatchHelper<vec_t, Vec<typename vec_t::scalar_t, 2>> {
157  typedef Vec<typename vec_t::scalar_t, 2> param_vec_t;
158  typedef typename vec_t::scalar_t scalar_t;
159  enum {dim = vec_t::dim,
160  proj_dim = dim + 1,
161  param_dim = param_vec_t::dim};
162  typedef Vec<scalar_t, proj_dim> proj_vec_t;
163 
164  static Range<Range<proj_vec_t>> cpToWcp(const Range<Range<vec_t>>& cp,
165  const Range<Range<double>>& w) {
166  Range<Range<proj_vec_t>> weighted_control_points;
167 
168  weighted_control_points.resize(cp.size());
169  for (int i = 0; i < cp.size(); i++) {
170  weighted_control_points[i].resize(cp[i].size());
171  for (int j = 0; j < cp[i].size(); j++) {
172  proj_vec_t temp;
173 
174  for (int y = 0; y < dim; y++) {
175  temp(y) = w[i][j] * cp[i][j](y);
176  }
177  temp(proj_dim - 1) = w[i][j];
178 
179  weighted_control_points[i][j] = temp;
180  }
181  }
182 
183  return weighted_control_points;
184  }
185 
186  static void computeDerivativeStructure(NURBSPatch<vec_t, param_vec_t>& patch) {
187  if (patch.der_structure != nullptr) {
188  return;
189  }
190 
191  // Calculate control points for the second dimension.
192  Range<Range<proj_vec_t>> v_der_wcp(patch.weighted_control_points.size());
193  for (int i = 0; i < patch.weighted_control_points.size(); i++) {
195  patch.weighted_control_points[i], patch.knots[1], v_der_wcp[i]);
196  }
197 
198  // Calculate knots for the second dimension.
199  Range<scalar_t> v_der_knots;
200  cad_helpers::generate_b_spline_derivative_knots(patch.knots[1], v_der_knots);
201 
202  // Calculate control points for the first dimension.
203  // Not efficient, but good enough, since the number of control points is usually
204  // sufficiently small. Optimization would require an additional data structure.
205  Range<Range<proj_vec_t>> u_der_wcp_initial(patch.weighted_control_points[0].size());
206  for (int i = 0; i < patch.weighted_control_points[0].size(); i++) {
207  Range<proj_vec_t> temp_range(patch.weighted_control_points.size());
208  for (int j = 0; j < patch.weighted_control_points.size(); j++) {
209  temp_range[j] = patch.weighted_control_points[j][i];
210  }
212  patch.knots[0], u_der_wcp_initial[i]);
213  }
214 
215  Range<Range<proj_vec_t>> u_der_wcp(patch.weighted_control_points.size() - 1);
216  for (int i = 0; i < patch.weighted_control_points.size() - 1; i++) {
217  u_der_wcp[i].resize(patch.weighted_control_points[0].size());
218  for (int j = 0; j < patch.weighted_control_points[0].size(); j++) {
219  u_der_wcp[i][j] = u_der_wcp_initial[j][i];
220  }
221  }
222 
223  // Calculate knots for the first dimension.
224  Range<scalar_t> u_der_knots;
225  cad_helpers::generate_b_spline_derivative_knots(patch.knots[0], u_der_knots);
226 
227  // Construct pair of derivative surfaces.
228  patch.der_structure =
229  std::unique_ptr<std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>>(
230  new std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>{
231  NURBSPatch<vec_t, param_vec_t>(u_der_wcp,
232  std::array<Range<scalar_t>, param_dim>{u_der_knots, patch.knots[1]},
233  std::array<int, param_dim>{patch.p[0] - 1, patch.p[1]}),
234  NURBSPatch<vec_t, param_vec_t>(v_der_wcp,
235  std::array<Range<scalar_t>, param_dim>{patch.knots[0], v_der_knots},
236  std::array<int, param_dim>{patch.p[0], patch.p[1] - 1})});
237  }
238 
239  static proj_vec_t evaluateWeighted(const NURBSPatch<vec_t, param_vec_t>& patch,
240  const param_vec_t& t) {
241  // This can maybe be calculated faster. We are especially limited by evaluating in
242  // arbitrary points, since NURBS evaluation can be done much faster on a grid.
243 
244  // Evaluate along the second dimension.
245  Range<proj_vec_t> secondary_control_points(patch.weighted_control_points.size());
246  for (int i = 0; i < patch.weighted_control_points.size(); i++) {
247  secondary_control_points[i] = cad_helpers::evaluate_b_spline(t(1), patch.p[1],
248  patch.weighted_control_points[i], patch.knots[1]);
249  }
250 
251  // Evaluate along the first dimension.
252  return cad_helpers::evaluate_b_spline(t(0), patch.p[0], secondary_control_points,
253  patch.knots[0]);
254  }
255 
256  static Eigen::Matrix<scalar_t, dim, param_dim> jacobian(
257  const NURBSPatch<vec_t, param_vec_t>& patch, const param_vec_t& t, const vec_t& pt,
258  const scalar_t& w) {
259  assert_msg(patch.der_structure != nullptr, "Derivative structure not initialised, "
260  "call initDerivative() before calling this function.");
261 
262  // Evaluate point on the first derivative NURBS surface.
263  proj_vec_t der_u_pt = (*patch.der_structure)[0].evaluateWeighted(t);
264 
265  // Calculate first partial derivative.
266  vec_t der_u;
267  for (int i = 0; i < dim; i++) {
268  der_u(i) = (der_u_pt(i) - der_u_pt(proj_dim - 1) * pt(i)) / w;
269  }
270 
271  // Evaluate point on the first derivative NURBS surface.
272  proj_vec_t der_v_pt = (*patch.der_structure)[1].evaluateWeighted(t);
273 
274  // Calculate second partial derivative.
275  vec_t der_v;
276  for (int i = 0; i < dim; i++) {
277  der_v(i) = (der_v_pt(i) - der_v_pt(proj_dim - 1) * pt(i)) / w;
278  }
279 
280  // Construct Jacobian matrix.
281  Eigen::Matrix<scalar_t, dim, 2> jm;
282 
283  jm.col(0) << der_u;
284  jm.col(1) << der_v;
285 
286  return jm;
287  }
288 
289  static Range<NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>>> getBoundaries(
290  const NURBSPatch<vec_t, param_vec_t>& patch) {
291  // Boundaries along the v direction.
292  NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>> b1(patch.weighted_control_points.front(),
293  {patch.knots[1]}, {patch.p[1]});
294  NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>> b2(patch.weighted_control_points.back(),
295  {patch.knots[1]}, {patch.p[1]});
296 
297  // Recalculate the range of control to use along the first dimension.
298  Range<proj_vec_t> wcp_front(patch.weighted_control_points.size());
299  Range<proj_vec_t> wcp_back(patch.weighted_control_points.size());
300  for (int i = 0; i < patch.weighted_control_points.size(); i++) {
301  wcp_front[i] = patch.weighted_control_points[i].front();
302  wcp_back[i] = patch.weighted_control_points[i].back();
303  }
304 
305  // Boundaries along the u direction.
306  NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>> b3(wcp_front, {patch.knots[0]},
307  {patch.p[0]});
308  NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>> b4(wcp_back, {patch.knots[0]},
309  {patch.p[0]});
310 
311  return Range<NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>>>{b1, b2, b3, b4};
312  }
313 
314  static param_vec_t getPatchParameterFromBoundaryParameter(
315  const NURBSPatch<vec_t, param_vec_t>& patch, const Vec<scalar_t, param_dim - 1>& t,
316  int i, scalar_t epsilon) {
317  scalar_t eps_1 = epsilon * (patch.knots[0].back() - patch.knots[0].front());
318  scalar_t eps_2 = epsilon * (patch.knots[1].back() - patch.knots[1].front());
319 
320  if (i == 0) return Vec<scalar_t, 2>(patch.knots[0].front() + eps_1, t(0));
321  else if (i == 1) return Vec<scalar_t, 2>(patch.knots[0].back() - eps_1, t(0));
322  else if (i == 2) return Vec<scalar_t, 2>(t(0), patch.knots[1].front() + eps_2);
323  else if (i == 3) return Vec<scalar_t, 2>(t(0), patch.knots[1].back() - eps_2);
324 
325  assert_msg(false, "Invalid boundary index, should be 0 <= i <= 3, got %d.", i);
326  return {};
327  }
328 
329  static Vec<scalar_t, param_dim - 1> getBoundaryParameterFromPatchParameter(
330  const NURBSPatch<vec_t, param_vec_t>& patch, const param_vec_t& t, int i) {
331  if (i == 0 || i == 1) return Vec<scalar_t, param_dim - 1>(t(1));
332  else if (i == 2 || i == 3) return Vec<scalar_t, param_dim - 1>(t(0));
333 
334  assert_msg(false, "Invalid edge index, should be 0 <= i <= 3, got %d.", i);
335  return {};
336  }
337 };
338 
339 template <typename vec_t>
340 struct nurbs_patch_internal::NURBSPatchHelper<vec_t, Vec<typename vec_t::scalar_t, 1>> {
341  typedef Vec<typename vec_t::scalar_t, 1> param_vec_t;
342  typedef typename vec_t::scalar_t scalar_t;
343  enum {dim = vec_t::dim,
344  proj_dim = dim + 1,
345  param_dim = param_vec_t::dim};
346  typedef Vec<scalar_t, proj_dim> proj_vec_t;
347 
348  static Range<proj_vec_t> cpToWcp(const Range<vec_t>& cp, const Range<double>& w) {
349  Range<proj_vec_t> weighted_control_points;
350 
351  weighted_control_points.resize(cp.size());
352  for (int i = 0; i < cp.size(); i++) {
353  proj_vec_t temp;
354 
355  for (int j = 0; j < dim; j++) {
356  temp(j) = w[i] * cp[i](j);
357  }
358  temp(proj_dim - 1) = w[i];
359 
360  weighted_control_points[i] = temp;
361  }
362 
363  return weighted_control_points;
364  }
365 
366  static void computeDerivativeStructure(NURBSPatch<vec_t, param_vec_t>& patch) {
367  if (patch.der_structure != nullptr || patch.p[0] == 0) {
368  return;
369  }
370 
371  Range<proj_vec_t> derivative_weighted_control_points;
372  Range<scalar_t> derivative_knots;
373  cad_helpers::generate_b_spline_derivative(patch.p[0], patch.weighted_control_points,
374  patch.knots[0], derivative_weighted_control_points, derivative_knots);
375 
376  patch.der_structure =
377  std::unique_ptr<std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>>(
378  new std::array<NURBSPatch<vec_t, param_vec_t>, param_dim>{
379  NURBSPatch<vec_t, param_vec_t>(derivative_weighted_control_points,
380  {derivative_knots}, {patch.p[0] - 1})});
381  }
382 
383  static proj_vec_t evaluateWeighted(const NURBSPatch<vec_t, param_vec_t>& patch,
384  const param_vec_t& t) {
385  return cad_helpers::evaluate_b_spline(t(0), patch.p[0], patch.weighted_control_points,
386  patch.knots[0]);
387  }
388 
389  static Eigen::Matrix<scalar_t, dim, param_dim> jacobian(
390  const NURBSPatch<vec_t, param_vec_t>& patch, const param_vec_t& t, const vec_t& pt,
391  const scalar_t& w) {
392  if (patch.p[0] == 0) {
393  return vec_t();
394  }
395  assert_msg(patch.der_structure != nullptr, "Derivative structure not initialised, "
396  "call initDerivative() before calling this function.");
397 
398  // Evaluate point on the derivative NURBS.
399  proj_vec_t der_pt = (*patch.der_structure)[0].evaluateWeighted(t);
400 
401  // Calculate Jacobian
402  vec_t jacobian;
403  for (int i = 0; i < dim; i++) {
404  jacobian(i) = (der_pt(i) - der_pt(proj_dim - 1) * pt(i)) / w;
405  }
406 
407  return jacobian;
408  }
409 
410  static Range<NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>>> getBoundaries(
411  const NURBSPatch<vec_t, param_vec_t>& patch) {
412  return {NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>>({}, {}, {}),
413  NURBSPatch<vec_t, Vec<scalar_t, param_dim - 1>>({}, {}, {})};
414  }
415 
416  static param_vec_t getPatchParameterFromBoundaryParameter(
417  const NURBSPatch<vec_t, param_vec_t>& patch, const Vec<scalar_t, param_dim - 1>& t,
418  int i, scalar_t epsilon) {
419  BoxShape<param_vec_t> dom = patch.getDomain();
420  scalar_t eps = epsilon * (patch.knots[0].back() - patch.knots[0].front());
421 
422  if (i == 0) return dom.beg() + param_vec_t(eps);
423  else if (i == 1) return dom.end() - param_vec_t(eps);
424 
425  assert_msg(false, "Invalid boundary index, should be 0 <= i <= 1, got %d.", i);
426  return {};
427  }
428 
429  static Vec<scalar_t, param_dim - 1> getBoundaryParameterFromPatchParameter(
430  const NURBSPatch<vec_t, param_vec_t>& patch, const param_vec_t& t, int i) {
431  if (i == 0 || i == 1) return {};
432 
433  assert_msg(false, "Invalid edge index, should be 0 <= i <= 1, got %d.", i);
434  return {};
435  }
436 };
438 
439 } // namespace mm
440 
441 #endif // MEDUSA_BITS_DOMAINS_NURBSPATCH_HPP_
mm::NURBSPatch::getPatchParameterFromBoundaryParameter
param_vec_t getPatchParameterFromBoundaryParameter(const Vec< scalar_t, param_dim - 1 > &t, int i, scalar_t epsilon=0) const
Get vector from the parametric domain of the NURBS patch from the parameter of the NURBS patch repres...
Definition: NURBSPatch.hpp:140
mm
Root namespace for the whole library.
Definition: Gaussian.hpp:14
scalar_t
Scalar scalar_t
Type of the elements, alias of Scalar.
Definition: MatrixBaseAddons.hpp:16
mm::Vec
Eigen::Matrix< scalar_t, dim, 1, Eigen::ColMajor|Eigen::AutoAlign, dim, 1 > Vec
Fixed size vector type, representing a mathematical 1d/2d/3d vector.
Definition: Vec_fwd.hpp:31
mm::cad_helpers::generate_b_spline_derivative_control_points
void generate_b_spline_derivative_control_points(int p, const Range< Vec< scalar_t, dim >> &control_points, const Range< scalar_t > &knots, Range< Vec< scalar_t, dim >> &der_control_points)
Generate control points of a B-spline that is the first derivative of the inputed B-spline.
Definition: cad_helpers.hpp:69
mm::NURBSPatch::evaluateWeighted
proj_vec_t evaluateWeighted(const param_vec_t &t) const
Evaluate NURBS patch in one point.
Definition: NURBSPatch.hpp:88
mm::NURBSPatch::NURBSPatch
NURBSPatch(const NdRange< proj_vec_t > &wcp, const std::array< Range< scalar_t >, param_dim > &ks, const std::array< int, param_dim > &in_p)
Construct NURBS patch with weighted control points.
Definition: NURBSPatch.hpp:17
mm::NURBSPatch::evaluate
vec_t evaluate(const param_vec_t &t, scalar_t *w) const
Evaluate NURBS patch in one point along the second dimension.
Definition: NURBSPatch.hpp:69
mm::NURBSPatch::jacobian
Eigen::Matrix< scalar_t, dim, param_dim > jacobian(const param_vec_t &t, const vec_t &pt, const scalar_t &w) const
Evaluate NURBS jacobian matrix in one point.
Definition: NURBSPatch.hpp:95
dim
@ dim
Number of elements of this matrix.
Definition: MatrixBaseAddons.hpp:14
mm::nurbs_patch_internal::NURBSPatchHelper
Internal structure of NURBSPatch that helps with partial class specialization.
Definition: NURBSPatch_fwd.hpp:34
assert_msg
#define assert_msg(cond,...)
Assert with better error reporting.
Definition: assert.hpp:75
mm::cad_helpers::generate_b_spline_derivative
void generate_b_spline_derivative(int p, const Range< Vec< scalar_t, dim >> &control_points, const Range< scalar_t > &knots, Range< Vec< scalar_t, dim >> &der_control_points, Range< scalar_t > &der_knots)
Generate control points and knot vector of a B-spline that is the first derivative of the inputed B-s...
Definition: cad_helpers.hpp:61
mm::NURBSPatch::getBoundaryParameterFromPatchParameter
Vec< scalar_t, param_dim - 1 > getBoundaryParameterFromPatchParameter(const param_vec_t &t, int i) const
Get parameter from the parametric domain of the NURBS patch representing the boundary of the NURBS pa...
Definition: NURBSPatch.hpp:149
mm::NURBSPatch
Class representing a single NURBS patch in an arbitrary dimensional space, defined on an arbitrary di...
Definition: NURBSPatch_fwd.hpp:72
mm::NURBSPatch::scalar_t
vec_t::scalar_t scalar_t
Scalar type.
Definition: NURBSPatch_fwd.hpp:76
mm::NURBSPatch::computeDerivativeStructure
void computeDerivativeStructure()
Calculate data needed for derivative evaluation.
Definition: NURBSPatch.hpp:64
mm::NURBSPatch::getBoundaries
Range< NURBSPatch< vec_t, Vec< scalar_t, param_dim - 1 > > > getBoundaries() const
Get all boundaries of the NURBS patch as a Range of NURBS patches.
Definition: NURBSPatch.hpp:135
mm::NURBSPatch::getDomain
BoxShape< param_vec_t > getDomain() const
Get domain of the NURBS patch.
Definition: NURBSPatch.hpp:122
mm::cad_helpers::generate_b_spline_derivative_knots
void generate_b_spline_derivative_knots(const Range< scalar_t > &knots, Range< scalar_t > &der_knots)
Generate knots of a B-spline that is the first derivative of the inputed B-spline.
Definition: cad_helpers.hpp:85
mm::BoxShape
Class for working with box shaped domains.
Definition: BoxShape_fwd.hpp:28
mm::NURBSPatch::evaluatePointAndJacobian
std::pair< vec_t, Eigen::Matrix< scalar_t, dim, param_dim > > evaluatePointAndJacobian(const param_vec_t &t) const
Evaluate NURBS and its jacobian matrix in one parameter.
Definition: NURBSPatch.hpp:114
assert.hpp
NURBSPatch_fwd.hpp
mm::NURBSPatch::operator=
NURBSPatch< vec_t, param_vec_t > & operator=(const NURBSPatch< vec_t, param_vec_t > &) noexcept
Copy assignment.
Definition: NURBSPatch.hpp:47
mm::cad_helpers::evaluate_b_spline
Vec< scalar_t, dim > evaluate_b_spline(scalar_t t, int p, const Range< Vec< scalar_t, dim >> &control_points, const Range< scalar_t > &knots, int k)
Evaluate B-spline in one point using De Boor's algorithm - , where is the number of dimensions.
Definition: cad_helpers.hpp:18
mm::NURBSPatch::proj_vec_t
Vec< scalar_t, proj_dim > proj_vec_t
Vector type before projection.
Definition: NURBSPatch_fwd.hpp:81
mm::NURBSPatch::NdRange
typename nurbs_patch_internal::NestedRange< param_dim, T >::type NdRange
Range of param_dim dimensions.
Definition: NURBSPatch_fwd.hpp:85
cad_helpers.hpp
Matrix
Matrix(const Scalar &s)
Construct matrix from scalar. Enabled only for fixed size matrices.
Definition: MatrixAddons.hpp:21
mm::Range
An extension of std::vector<T> to support additional useful operations.
Definition: Range_fwd.hpp:30