NurbsUtils Namespace Reference

Classes

class  array2
 

Functions

bool isKnotVectorMonotonic (const std::vector< double > &knots)
 
bool close (double a, double b, double eps)
 
int findSpan (int degree, const std::vector< double > &knots, double u)
 
double bsplineOneBasis (int i, int deg, const std::vector< double > &U, double u)
 
void bsplineBasis (int deg, int span, const std::vector< double > &knots, double u, std::vector< double > &N)
 
void bsplineDerBasis (int deg, int span, const std::vector< double > &knots, double u, int nDers, std::vector< std::vector< double >> &ders)
 
std::vector< MdoublecreateUniformKnotVector (unsigned int numberOfControlPoints, unsigned int degree, bool clampedAtStart, bool clampedAtEnd)
 Creates a uniform (clamped) knot vector. More...
 
void normalizeKnotVector (std::vector< Mdouble > &knots)
 Resets the knot vector to the interval [0, 1]. More...
 
void extendKnotVector (std::vector< Mdouble > &knots, unsigned int degree, unsigned int numStart, unsigned int numEnd, bool forceBothEndsUniform=false)
 Extends the knot vector for when control points have been added at the start or end. More...
 
Vec3D evaluate (Mdouble u, Mdouble v, std::vector< Mdouble > knotsU, std::vector< Mdouble > knotsV, std::vector< std::vector< Vec3D >> controlPoints, std::vector< std::vector< Mdouble >> weights)
 Evaluate point on a NURBS surface. More...
 

Function Documentation

◆ bsplineBasis()

void NurbsUtils::bsplineBasis ( int  deg,
int  span,
const std::vector< double > &  knots,
double  u,
std::vector< double > &  N 
)

Compute all non-zero B-spline basis functions

Parameters
[in]degDegree of the basis function.
[in]spanIndex obtained from findSpan() corresponding the u and knots.
[in]knotsKnot vector corresponding to the basis functions.
[in]uParameter to evaluate the basis functions at.
[in,out]NValues of (deg+1) non-zero basis functions.
131 {
132  N.clear();
133  N.resize(deg + 1, 0.0);
134  std::vector<double> left, right;
135  left.resize(deg + 1, 0.0);
136  right.resize(deg + 1, 0.0);
137 
138  N[0] = 1.0;
139 
140  for (int j = 1; j <= deg; j++)
141  {
142  left[j] = (u - knots[span + 1 - j]);
143  right[j] = knots[span + j] - u;
144  Mdouble saved = 0.0;
145  for (int r = 0; r < j; r++)
146  {
147  const Mdouble temp = N[r] / (right[r + 1] + left[j - r]);
148  N[r] = saved + right[r + 1] * temp;
149  saved = left[j - r] * temp;
150  }
151  N[j] = saved;
152  }
153 }
double Mdouble
Definition: GeneralDefine.h:34

Referenced by NurbsSurface::evaluate(), and evaluate().

◆ bsplineDerBasis()

void NurbsUtils::bsplineDerBasis ( int  deg,
int  span,
const std::vector< double > &  knots,
double  u,
int  nDers,
std::vector< std::vector< double >> &  ders 
)

Compute all non-zero derivatives of B-spline basis functions

Parameters
[in]degDegree of the basis function.
[in]spanIndex obtained from findSpan() corresponding the u and knots.
[in]knotsKnot vector corresponding to the basis functions.
[in]uParameter to evaluate the basis functions at.
[in]nDersNumber of derivatives to compute (nDers <= deg)
[in,out]dersValues of non-zero derivatives of basis functions.
156  {
157 
158  std::vector<double> left, right;
159  left.resize(deg + 1, 0.0);
160  right.resize(deg + 1, 0.0);
161  double saved = 0.0, temp = 0.0;
162 
163  array2<double> ndu(deg + 1, deg + 1);
164  ndu(0, 0) = 1.0;
165 
166  for (int j = 1; j <= deg; j++) {
167  left[j] = u - knots[span + 1 - j];
168  right[j] = knots[span + j] - u;
169  saved = 0.0;
170 
171  for (int r = 0; r < j; r++) {
172  // Lower triangle
173  ndu(j, r) = right[r + 1] + left[j - r];
174  temp = ndu(r, j - 1) / ndu(j, r);
175  // Upper triangle
176  ndu(r, j) = saved + right[r + 1] * temp;
177  saved = left[j - r] * temp;
178  }
179 
180  ndu(j, j) = saved;
181  }
182 
183  ders.clear();
184  ders.resize(nDers + 1);
185  for (int i = 0; i < ders.size(); i++) {
186  ders[i].resize(deg + 1, 0.0);
187  }
188 
189  for (int j = 0; j <= deg; j++) {
190  ders[0][j] = ndu(j, deg);
191  }
192 
193  array2<double> a(2, deg + 1);
194 
195  for (int r = 0; r <= deg; r++) {
196  int s1 = 0;
197  int s2 = 1;
198  a(0, 0) = 1.0;
199 
200  for (int k = 1; k <= nDers; k++) {
201  double d = 0.0;
202  int rk = r - k;
203  int pk = deg - k;
204  int j1 = 0;
205  int j2 = 0;
206 
207  if (r >= k) {
208  a(s2, 0) = a(s1, 0) / ndu(pk + 1, rk);
209  d = a(s2, 0) * ndu(rk, pk);
210  }
211 
212  if (rk >= -1) {
213  j1 = 1;
214  }
215  else {
216  j1 = -rk;
217  }
218 
219  if (r - 1 <= pk) {
220  j2 = k - 1;
221  }
222  else {
223  j2 = deg - r;
224  }
225 
226  for (int j = j1; j <= j2; j++) {
227  a(s2, j) = (a(s1, j) - a(s1, j - 1)) / ndu(pk + 1, rk + j);
228  d += a(s2, j) * ndu(rk + j, pk);
229  }
230 
231  if (r <= pk) {
232  a(s2, k) = -a(s1, k - 1) / ndu(pk + 1, r);
233  d += a(s2, k) * ndu(r, pk);
234  }
235 
236 
237  ders[k][r] = d;
238 
239  int temp = s1;
240  s1 = s2;
241  s2 = temp;
242  }
243  }
244 
245  double fac = static_cast<double>(deg);
246  for (int k = 1; k <= nDers; k++) {
247  for (int j = 0; j <= deg; j++) {
248  ders[k][j] *= fac;
249  }
250  fac *= static_cast<double>(deg - k);
251  }
252 }
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References constants::i.

Referenced by NurbsSurface::evaluateDerivatives().

◆ bsplineOneBasis()

double NurbsUtils::bsplineOneBasis ( int  i,
int  deg,
const std::vector< double > &  U,
double  u 
)

Compute a single B-spline basis function

Parameters
[in]iThe ith basis function to compute.
[in]degDegree of the basis function.
[in]UKnot vector corresponding to the basis functions.
[in]uParameter to evaluate the basis functions at.
Returns
The value of the ith basis function at u.
85 {
86  int m = static_cast<int>(U.size()) - 1;
87  // Special case
88  if ((i == 0 && close(u, U[0])) || (i == m - deg - 1 && close(u, U[m])))
89  {
90  return 1.0;
91  }
92  // Local Property
93  if (u < U[i] || u >= U[i + deg + 1])
94  {
95  return 0.0;
96  }
97  // Initialize zeroth-degree functions
98  std::vector<double> N;
99  N.resize(deg + 1);
100  for (int j = 0; j <= deg; j++)
101  {
102  N[j] = (u >= U[i + j] && u < U[i + j + 1]) ? 1.0 : 0.0;
103  }
104  // Compute triangular table
105  for (int k = 1; k <= deg; k++)
106  {
107  double saved = (close(N[0], 0.0)) ? 0.0
108  : ((u - U[i]) * N[0]) / (U[i + k] - U[i]);
109  for (int j = 0; j < deg - k + 1; j++)
110  {
111  double Uleft = U[i + j + 1];
112  double Uright = U[i + j + k + 1];
113  if (close(N[j + 1], 0.0))
114  {
115  N[j] = saved;
116  saved = 0.0;
117  }
118  else
119  {
120  double temp = N[j + 1] / (Uright - Uleft);
121  N[j] = saved + (Uright - u) * temp;
122  saved = (u - Uleft) * temp;
123  }
124  }
125  }
126  return N[0];
127 }
bool close(double a, double b, double eps)
Definition: NurbsUtils.cc:38

References close(), and constants::i.

◆ close()

bool NurbsUtils::close ( double  a,
double  b,
double  eps 
)
39 {
40  return (std::abs(a - b) < eps) ? true : false;
41 }

Referenced by bsplineOneBasis(), extendKnotVector(), findSpan(), and normalizeKnotVector().

◆ createUniformKnotVector()

std::vector< Mdouble > NurbsUtils::createUniformKnotVector ( unsigned int  numberOfControlPoints,
unsigned int  degree,
bool  clampedAtStart,
bool  clampedAtEnd 
)

Creates a uniform (clamped) knot vector.

Parameters
numberOfControlPointsThe number of control points the knot vector should correspond to
degreeThe degree the knot vector should correspond to
clampedAtStartWhether or not the knot vector should be clamped at the start
clampedAtEndWhether or not the knot vector should be clamped at the end
Returns
A uniform (clamped) knot vector
255 {
256  // Total number of knots
257  unsigned int numberOfKnots = numberOfControlPoints + degree + 1;
258  std::vector<Mdouble> knots;
259  knots.reserve(numberOfKnots);
260 
261  // Create uniform from 0 to 1
262  for (int i = 0; i < numberOfKnots; i++)
263  {
264  knots.push_back(i / static_cast<Mdouble>(numberOfKnots - 1));
265  }
266 
267  // When clamped at start, first degree+1 values are set to 0.
268  if (clampedAtStart)
269  {
270  for (int i = 0; i <= degree; i++)
271  {
272  knots[i] = 0.0;
273  }
274  }
275 
276  // When clamped at end, last degree+1 values are set to 1
277  if (clampedAtEnd)
278  {
279  for (int i = 0; i <= degree; i++)
280  {
281  knots[knots.size() - 1 - i] = 1.0;
282  }
283  }
284 
285  return knots;
286 }
const Mdouble degree
Definition: ExtendedMath.h:52

References constants::degree, and constants::i.

Referenced by NurbsSurface::NurbsSurface().

◆ evaluate()

Vec3D NurbsUtils::evaluate ( Mdouble  u,
Mdouble  v,
std::vector< Mdouble knotsU,
std::vector< Mdouble knotsV,
std::vector< std::vector< Vec3D >>  controlPoints,
std::vector< std::vector< Mdouble >>  weights 
)

Evaluate point on a NURBS surface.

Parameters
uParameter to evaluate the surface at
vParameter to evaluate the surface at
knotsUKnot vector in u-direction
knotsVKnot vector in v-direction
controlPoints2D vector of control points
weights2D vector of weights
Returns
Resulting point on the surface at (u, v)
378 {
379  unsigned int degreeU = knotsU.size() - controlPoints.size() - 1;
380  unsigned int degreeV = knotsV.size() - controlPoints[0].size() - 1;
381 
382  Vec3D point = {0,0,0};
383  double temp = 0;
384 
385  // Find span and non-zero basis functions
386  int spanU = findSpan(degreeU, knotsU, u);
387  int spanV = findSpan(degreeV, knotsV, v);
388  std::vector<double> Nu, Nv;
389  bsplineBasis(degreeU, spanU, knotsU, u, Nu);
390  bsplineBasis(degreeV, spanV, knotsV, v, Nv);
391  // make linear combination
392  for (int l = 0; l <= degreeV; ++l)
393  {
394  for (int k = 0; k <= degreeU; ++k)
395  {
396  double weight = Nv[l] * Nu[k] * weights[spanU - degreeU + k][spanV - degreeV + l];
397  point += weight * controlPoints[spanU - degreeU + k][spanV - degreeV + l];
398  temp += weight;
399  }
400  }
401  point /= temp;
402  return point;
403 }
Definition: Vector.h:51
double Nu
Poisson's ratio.
Definition: TwenteMeshGluing.cpp:67
int findSpan(int degree, const std::vector< double > &knots, double u)
Definition: NurbsUtils.cc:43
void bsplineBasis(int deg, int span, const std::vector< double > &knots, double u, std::vector< double > &N)
Definition: NurbsUtils.cc:129

References bsplineBasis(), findSpan(), and Global_Physical_Variables::Nu.

Referenced by WearableNurbsWall::getVolumeUnderSurface(), WearableNurbsWall::getVolumeUnderSurfaceX(), main(), and NurbsSurface::set().

◆ extendKnotVector()

void NurbsUtils::extendKnotVector ( std::vector< Mdouble > &  knots,
unsigned int  degree,
unsigned int  numStart,
unsigned int  numEnd,
bool  forceBothEndsUniform = false 
)

Extends the knot vector for when control points have been added at the start or end.

Parameters
knotsThe knot vector to be extended
degreeThe nurbs degree
numStartNumber to be added at the start
numEndNumber to be added at the end
forceBothEndsUniformWhen only adding to one end, the other end can remain untouched or it can be forced to be made uniform
305 {
306  // Usually it should already be in normalized form, but just to be sure.
307  // Needed because otherwise comparing to uniform step size doesn't work.
308  normalizeKnotVector(knots);
309 
310  Mdouble uniformStepSize = 1.0 / (knots.size() - 1);
311 
312  // This only works properly when the original start and end knots are uniformly spaced.
313  // Therefore force the degree+1 number of knots at the start and end to have uniform step size
314  // and remember if they weren't already uniform, so when needed a message can be logged.
315 
316  // Since the original knot vector might change a bit, only do stuff when there are actually knots to be added.
317  if (numStart > 0 || forceBothEndsUniform)
318  {
319  // Check if first degree+1 knots are uniform.
320  // Update: also an additional number of knots added to the end should be uniform (the ones "wrapping around")
321  bool startsOfUniform = true;
322  for (int i = 0; i <= degree + numEnd; i++)
323  {
324  if (!close(knots[i], i * uniformStepSize))
325  {
326  startsOfUniform = false;
327  knots[i] = i * uniformStepSize;
328  }
329  }
330 
331  // When needed, let user know that original knot vector has been edited.
332  if (!startsOfUniform)
333  {
334  logger(INFO, "Start of knot vector (first degree+1 values) has been changed to be uniform. The shape might be slightly affected. ");
335  }
336 
337  // Insert knots at start
338  for (int i = 1; i <= numStart; i++)
339  {
340  knots.insert(knots.begin(), -i * uniformStepSize);
341  }
342  }
343 
344  // Since the original knot vector might change a bit, only do stuff when there are actually knots to be added.
345  if (numEnd > 0 || forceBothEndsUniform)
346  {
347  // Check if last degree+1 knots are uniform.
348  // Update: also an additional number of knots added to the start should be uniform (the ones "wrapping around")
349  bool endsOfUniform = true;
350  for (int i = 0; i <= degree + numStart; i++)
351  {
352  if (!close(knots[knots.size() - 1 - i], (1.0 - i * uniformStepSize)))
353  {
354  endsOfUniform = false;
355  knots[knots.size() - 1 - i] = 1.0 - i * uniformStepSize;
356  }
357  }
358 
359  // When needed, let user know that original knot vector has been edited.
360  if (!endsOfUniform)
361  {
362  logger(INFO, "End of knot vector (last degree+1 values) has been changed to be uniform. The shape might be slightly affected. ");
363  }
364 
365  // Add knots to end
366  for (int i = 1; i <= numEnd; i++)
367  {
368  knots.push_back(1.0 + i * uniformStepSize);
369  }
370  }
371 
372  // Reset to interval [0, 1]
373  normalizeKnotVector(knots);
374 }
LL< Log::INFO > INFO
Info log level.
Definition: Logger.cc:55
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
void normalizeKnotVector(std::vector< Mdouble > &knots)
Resets the knot vector to the interval [0, 1].
Definition: NurbsUtils.cc:288

References close(), constants::degree, constants::i, INFO, logger, and normalizeKnotVector().

Referenced by NurbsSurface::wrapAroundInU(), and NurbsSurface::wrapAroundInV().

◆ findSpan()

int NurbsUtils::findSpan ( int  degree,
const std::vector< double > &  knots,
double  u 
)

Find the span of the given parameter in the knot vector.

Parameters
[in]degreeDegree of the curve.
[in]knotsKnot vector of the curve.
[in]uParameter value.
Returns
Span index into the knot vector such that (span - 1) < u <= span
44 {
45  // index of last control point
46  int n = static_cast<int>(knots.size()) - degree - 2;
47 
48  // For u that is equal to last knot value
49  if (close(u, knots[n + 1]))
50  {
51  return n;
52  }
53 
54  // For values of u that lies outside the domain
55  if (u > knots[n + 1])
56  {
57  return n;
58  }
59  if (u < knots[degree])
60  {
61  return degree;
62  }
63 
64  // Binary search
65  // TODO: Replace this with std::lower_bound
66  int low = degree;
67  int high = n + 1;
68  int mid = (int) std::floor((low + high) / 2.0);
69  while (u < knots[mid] || u >= knots[mid + 1])
70  {
71  if (u < knots[mid])
72  {
73  high = mid;
74  }
75  else
76  {
77  low = mid;
78  }
79  mid = (low + high) / 2;
80  }
81  return mid;
82 }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:32

References close(), constants::degree, and n.

Referenced by NurbsSurface::evaluate(), evaluate(), and NurbsSurface::evaluateDerivatives().

◆ isKnotVectorMonotonic()

bool NurbsUtils::isKnotVectorMonotonic ( const std::vector< double > &  knots)
34 {
35  return std::is_sorted(knots.begin(), knots.end());
36 }

Referenced by NurbsSurface::set().

◆ normalizeKnotVector()

void NurbsUtils::normalizeKnotVector ( std::vector< Mdouble > &  knots)

Resets the knot vector to the interval [0, 1].

Parameters
knotsThe knot vector to be normalized
289 {
290  // Reset the interval to [0,1]
291  const Mdouble minK = knots.front();
292  const Mdouble maxK = knots.back();
293 
294  // Already in normalized form
295  if (close(minK, 0.0) && close(maxK, 1.0))
296  return;
297 
298  for (Mdouble& k : knots)
299  {
300  k = (k - minK) / (maxK - minK);
301  }
302 }

References close().

Referenced by extendKnotVector(), NurbsSurface::set(), and NurbsSurface::unclampKnots().