NURBS3DCurve.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2007-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "NURBS3DCurve.H"
31 #include "vectorField.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 label NURBS3DCurve::sgn(const scalar val) const
41 {
42  return val>=scalar(0) ? 1 : -1;
43 }
44 
45 
46 scalar NURBS3DCurve::abs(const scalar val) const
47 {
48  return (sgn(val) == 1) ? val : -val;
49 }
50 
51 
52 label NURBS3DCurve::mod(const label x, const label interval) const
53 {
54  label ratio(x%interval);
55  return ratio < 0 ? ratio + interval : ratio;
56 }
57 
58 
59 void NURBS3DCurve::setUniformU()
60 {
61  const vectorField& curve(*this);
62  label nPts(curve.size());
63 
64  forAll(curve, ptI)
65  {
66  u_[ptI] = scalar(ptI)/scalar(nPts - 1);
67  }
68 }
69 
70 
71 bool NURBS3DCurve::bound
72 (
73  scalar& u,
74  const scalar minVal = 1e-7,
75  const scalar maxVal = 0.999999
76 ) const
77 {
78  // lower value bounding
79  if (u < scalar(0))
80  {
81  u = minVal;
82  return true;
83  //Info<< "Lower bound hit." << endl;
84  }
85 
86  // upper value bounding
87  if (u > scalar(1))
88  {
89  u = maxVal;
90  return true;
91  //Info<< "Upper bound hit." << endl;
92  }
93 
94  return false;
95 }
96 
97 
98 void NURBS3DCurve::setEquidistantU
99 (
100  scalarList& U,
101  const label lenAcc = 25,
102  const label maxIter = 10,
103  const label spacingCorrInterval=-1,
104  const scalar tolerance = 1.e-5
105 ) const
106 {
107  const label nPts(U.size());
108  const scalar xLength(length() /(nPts - 1));
109  const scalar uLength(scalar(1) / scalar(nPts - 1));
110 
111  U[0] = Zero;
112  U[nPts - 1] = scalar(1);
113 
114  for (label ptI = 1; ptI<(nPts - 1); ptI++)
115  {
116  const scalar UPrev(U[ptI - 1]);
117  scalar& UCurr(U[ptI]);
118  scalar direc(scalar(1));
119  scalar xDiff(scalar(0));
120  scalar delta(scalar(0));
121  bool overReached(false);
122 
123  UCurr = UPrev + uLength;
124 
125  // Find the starting U value to ensure target is within 1 uLength.
126  while (true)
127  {
128  bool bounded(bound(UCurr));
129 
130  delta = length(UPrev, UCurr, lenAcc);
131  xDiff = xLength - delta;
132 
133  // Found the point.
134  if (abs(xDiff) < tolerance)
135  {
136  break;
137  }
138  else
139  {
140  direc = sgn(xDiff);
141 
142  if (bounded && (direc == 1))
143  {
144  // uLength addition makes U exceed 1 so it becomes bounded.
145  // However, the desired x location still exceeds how far the
146  // bounded uLength can move (~e-5 error).
147  // Must force U to be u=0.999999.
148  overReached = true;
149  break;
150  }
151  else if (direc == scalar(1))
152  {
153  UCurr += uLength;
154  }
155  else
156  {
157  break;
158  }
159  }
160  }
161 
162  if (!overReached)
163  {
164  label iter(0);
165 
166  while (iter < maxIter)
167  {
168  // Set the new search length to ensure convergence and next v.
169  direc /= scalar(2);
170  UCurr += direc * uLength;
171 
172  bound(UCurr);
173 
174  // Can employ an occasional tolerance check from beg of curve.
175  if (
176  (spacingCorrInterval != -1)
177  && (mod(ptI, spacingCorrInterval) == 0)
178  )
179  {
180  delta = length(scalar(0), UCurr, ptI*lenAcc);
181  xDiff = (ptI * xLength) - delta;
182  }
183  else
184  {
185  delta = length(UPrev, UCurr, lenAcc);
186  xDiff = xLength - delta;
187  }
188 
189  // Break if found point or find the new search direction.
190  if (abs(xDiff) < tolerance)
191  {
192  break;
193  }
194  else
195  {
196  const scalar oldDirec(direc);
197  direc = sgn(xDiff) * abs(oldDirec);
198  }
199 
200  iter++;
201  }
202  }
203  }
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
208 
210 (
211  const NURBSbasis& basis,
212  const List<vector>& CPs,
213  const List<scalar>& weights,
214  const scalarField& u,
215  const label nPts,
216  const word name
217 )
218 :
219  vectorField(nPts, Zero),
220  CPs_(CPs),
221  weights_(weights),
222  u_(u),
223  name_(name),
224  basis_(basis),
225 
226  givenInitNrm_(Zero),
227 
228  nrmOrientation_(ALIGNED)
229 
230 {
231  buildCurve();
232 }
233 
234 
236 (
237  const NURBSbasis& basis,
238  const List<vector>& CPs,
239  const List<scalar>& weights,
240  const label nPts,
241  const word name
242 )
243 :
244  vectorField(nPts, Zero),
245  CPs_(CPs),
246  weights_(weights),
247  u_(nPts, Zero),
248  name_(name),
249  basis_(basis),
250 
251  givenInitNrm_(Zero),
252 
253  nrmOrientation_(ALIGNED)
254 
255 {
256  setUniformU();
257  buildCurve();
258 }
259 
260 
262 (
263  const NURBSbasis& basis,
264  const List<vector>& CPs,
265  const label nPts,
266  const word name
267 )
268 :
269  vectorField(nPts, Zero),
270  CPs_(CPs),
271  weights_(CPs.size(), scalar(1)),
272  u_(nPts, Zero),
273  name_(name),
274  basis_(basis),
275 
276  givenInitNrm_(Zero),
277 
278  nrmOrientation_(1)
279 
280 {
281  setUniformU();
282  buildCurve();
283 }
284 
285 
286 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
287 
288 // * * * * * * * * * * * * * * * * Set Functions * * * * * * * * * * * * * * //
289 
291 (
292  const vector& givenNrm,
293  const vector& givenTan
294 )
295 {
296  givenInitNrm_ = givenNrm;
297 
298  const vector tan(curveDerivativeU(Zero));
299  vector curveNrm(tan ^ givenTan);
300 
301  if ((givenNrm & curveNrm) >= scalar(0))
302  {
303  nrmOrientation_ = ALIGNED;
304  }
305  else
306  {
307  nrmOrientation_ = OPPOSED;
308  }
309 
310  Info<< "Initial nrmOrientation after comparison to NURBS u = 0 nrm : "
311  << nrmOrientation_
312  << endl;
313 }
314 
315 
317 (
318  const vector& givenNrm,
319  const scalar zVal
320 )
321 {
322  givenInitNrm_ = givenNrm;
323 
324  const vector tan(curveDerivativeU(Zero));
325  vector curveNrm(Zero);
326 
327  curveNrm.x() = -tan.y();
328  curveNrm.y() = tan.x();
329  curveNrm.z() = zVal;
330 
331  if ((givenNrm & curveNrm) >= scalar(0))
332  {
333  nrmOrientation_ = ALIGNED;
334  }
335  else
336  {
337  nrmOrientation_ = OPPOSED;
338  }
339 
340  Info<< "Initial nrmOrientation after comparison to NURBS u = 0 nrm : "
341  << nrmOrientation_
342  << endl;
343 }
344 
345 
347 {
348  if (nrmOrientation_ == ALIGNED)
349  {
350  nrmOrientation_ = OPPOSED;
351  }
352  else
353  {
354  nrmOrientation_ = ALIGNED;
355  }
356 }
357 
358 
360 {
361  CPs_ = CPs;
362 }
363 
364 
366 {
367  weights_ = weights;
368 }
369 
370 
372 {
373  name_ = name;
374 }
375 
376 
378 {
379  const label degree(basis_.degree());
380 
381  // Iterate over the CPs of this curve.
382  forAll(*this, ptI)
383  {
384  this->operator[](ptI) = vector::zero;
385 
386  const scalar u(u_[ptI]);
387 
388  // Compute denominator.
389  scalar denom(Zero);
390 
391  forAll(CPs_, CPJ)
392  {
393  denom += basis_.basisValue(CPJ, degree, u)*weights_[CPJ];
394  }
395 
396  // Compute curve point.
397  forAll(CPs_, CPI)
398  {
399  this->operator[](ptI)
400  += CPs_[CPI]
401  * basis_.basisValue(CPI, degree, u)
402  * weights_[CPI]/denom;
403  }
404  }
405 }
406 
407 
409 {
410  Info<< "Inverting NURBS curve " << name_ << endl;
411 
412  const label nCPs(CPs_.size());
413  List<vector> invertedCPs(nCPs, Zero);
414  List<scalar> invertedWeights(nCPs, Zero);
415 
416  for (label CPI = 0; CPI<nCPs; CPI++)
417  {
418  invertedCPs[CPI] = CPs_[nCPs - 1 - CPI];
419  invertedWeights[CPI] = weights_[nCPs - 1 - CPI];
420  }
421 
422  CPs_ = invertedCPs;
423  weights_ = invertedWeights;
424 
425  buildCurve();
426 }
427 
428 
430 (
431  const label lenAcc,
432  const label maxIter,
433  const label spacingCorrInterval,
434  const scalar tolerance
435 )
436 {
437 /*
438  Info<< "Making points equidistant is physical space on curve "
439  << name_
440  << endl;
441 */
442  setEquidistantU
443  (
444  u_,
445  lenAcc,
446  maxIter,
447  spacingCorrInterval,
448  tolerance
449  );
450 
451  buildCurve();
452 }
453 
454 
455 // * * * * * * * * * * * * * Point Calc Functions * * * * * * * * * * * * * //
456 
457 vector NURBS3DCurve::curvePoint(const scalar u) const
458 {
459  // Compute denominator.
460  const label degree(basis_.degree());
461  scalar NW(Zero);
462 
463  forAll(CPs_, CPI)
464  {
465  NW += basis_.basisValue(CPI, degree, u)*weights_[CPI];
466  }
467 
468  // Compute curve point.
469  vector point(Zero);
470 
471  forAll(CPs_, CPI)
472  {
473  point +=
474  CPs_[CPI]
475  *basis_.basisValue(CPI, degree, u)
476  *weights_[CPI]/NW;
477  }
478 
479  return point;
480 }
481 
482 
484 (
485  const vector& targetPoint,
486  const label maxIter,
487  const scalar tolerance
488 )
489 {
490  // Loop through curve points to find the closest one for initialization.
491  const vectorField& curve(*this);
492  scalar dist(GREAT);
493  label closeI(-1);
494 
495  forAll(curve, ptI)
496  {
497  const scalar distLoc(mag(targetPoint - curve[ptI]));
498 
499  if (distLoc < dist)
500  {
501  dist = distLoc;
502  closeI = ptI;
503  }
504  }
505 
506  label iter(0);
507  scalar u(scalar(closeI)/scalar(this->size()-1));
508  vector xu(curvePoint(u));
509  scalar res(GREAT);
510 
511  do
512  {
513  vector dxdu(curveDerivativeU(u));
514  vector d2xdu2(curveDerivativeUU(u));
515  scalar lhs((dxdu&dxdu) + ((xu - targetPoint) & d2xdu2));
516  scalar rhs(-((xu - targetPoint) & dxdu));
517 
518  // Update parametric coordinate and compute new point and
519  // bound param coordinates if needed.
520  // Compute residual.
521  u += rhs/lhs;
522 
523  bound(u);
524 
525  xu = curvePoint(u);
526  res = mag((xu - targetPoint) & curveDerivativeU(u));
527 
528  //Info<< "dxdu" << dxdu << endl;
529  //Info<< "d2xdu2" << d2xdu2 << endl;
530  //Info<< "Iteration " << iter << ", Residual " << res << endl;
531  }
532  while ((iter++< maxIter) && (res > tolerance));
533 /*
534  // debug info
535  Info<< "Residual after " << iter << " iterations : " << res
536  << "\nparametric coordinate " << u
537  << "\ncurve point closest to target " << xu
538  << "\ntarget being " << targetPoint
539  << endl;
540 */
541  // warning if method did not reach threshold
542  if (iter > maxIter)
543  {
545  << "Finding curve point closest to " << targetPoint << " failed."
546  << endl;
547  }
548 
549  return u;
550 }
551 
552 
554 (
555  const vector& targetPoint,
556  const scalar initGuess,
557  const label maxIter,
558  const scalar tolerance
559 )
560 {
561  // Loop through curve points to find the closest one for initialization.
562  label iter(0);
563  scalar u(initGuess);
564  vector xu(curvePoint(u));
565  scalar res(GREAT);
566 
567  do
568  {
569  vector dxdu(curveDerivativeU(u));
570  vector d2xdu2(curveDerivativeUU(u));
571  scalar lhs((dxdu&dxdu) + ((xu - targetPoint) & d2xdu2));
572  scalar rhs(-((xu - targetPoint) & dxdu));
573 
574  // Update parametric coordinate and compute new point and
575  // bound param coordinates if needed.
576  // Compute residual.
577  u += rhs/lhs;
578 
579  bound(u);
580 
581  xu = curvePoint(u);
582  res = mag((xu - targetPoint) & curveDerivativeU(u));
583 
584  //Info<< "dxdu" << dxdu << endl;
585  //Info<< "d2xdu2" << d2xdu2 << endl;
586  //Info<< "bound u " << u << endl;
587  //Info<< "u " << u << endl;
588  //Info<< "Iteration " << iter << ", Residual " << res << endl;
589  }
590  while ((iter++< maxIter) && (res > tolerance));
591 /*
592  // debug info
593  Info<< "Residual after " << iter << " iterations : " << res
594  << "\nparametric coordinate " << u
595  << "\ncurve point closest to target " << xu
596  << "\ntarget being " << targetPoint
597  << endl;
598 */
599  // warning if method did not reach threshold
600  if (iter > maxIter)
601  {
603  << "Finding curve point closest to " << targetPoint
604  << " failed."
605  << endl;
606  }
607 
608  return u;
609 }
610 
611 
612 const vector NURBS3DCurve::nrm3D(const vector& refTan, const scalar u) const
613 {
614  vector curveNrm(Zero);
615 
616  if (nrmOrientation_ == ALIGNED)
617  {
618  curveNrm = curveDerivativeU(u) ^ refTan;
619  }
620  else
621  {
622  curveNrm = refTan ^ curveDerivativeU(u);
623  }
624 
625  curveNrm.normalise();
626 
627  return curveNrm;
628 }
629 
630 
631 const vector NURBS3DCurve::nrm2D(const scalar zVal, const scalar u) const
632 {
633  const vector tan(curveDerivativeU(u));
634  vector curveNrm(Zero);
635 
636  curveNrm.x() = -nrmOrientation_*tan.y();
637  curveNrm.y() = nrmOrientation_*tan.x();
638  curveNrm.z() = zVal;
639  curveNrm /= mag(curveNrm);
640 
641  return curveNrm;
642 }
643 
644 
646 (
647  const scalarField& oldKnots,
648  const scalar uBar,
649  const label kInsert
650 )
651 {
652  // Get the req ref info.
653  // Insertion into curve of non-uniform weight is not currently supported.
654  const label degree(basis_.degree());
655  const label nCPs(basis_.nCPs());
656  List<vector> newCPs(nCPs, Zero);
657  List<scalar> newWeights(nCPs, scalar(1));
658 
659  // Compute the new CPs and their weights.
660  for (label CPI = 0; CPI < (kInsert - degree + 1); CPI++)
661  {
662  newCPs[CPI] = CPs_[CPI];
663  }
664 
665  for (label CPI = (kInsert - degree + 1); CPI < (kInsert + 1); CPI++)
666  {
667  const scalar uIOld(oldKnots[CPI]);
668  const scalar uIDOld(oldKnots[CPI + degree]);
669  const scalar ratio((uBar - uIOld) /(uIDOld - uIOld));
670 
671  newCPs[CPI] = (ratio*CPs_[CPI] + (1 - ratio)*CPs_[CPI - 1]);
672  }
673 
674  for (label CPI= (kInsert + 1); CPI<newCPs.size(); CPI++)
675  {
676  newCPs[CPI] = CPs_[CPI - 1];
677  }
678 
679  // Reset the CPs and weights and recompute the curve.
680  CPs_ = newCPs;
681  weights_ = newWeights;
682 
683  buildCurve();
684 }
685 
686 
688 (
689  const label nPts,
690  const label lenAcc,
691  const label maxIter,
692  const label spacingCorrInterval,
693  const scalar tolerance
694 )
695 {
696 /*
697  Info<< "Generating points equidistant in physical space on curve "
698  << name_
699  << endl;
700 */
701  scalarList U(nPts, Zero);
702 
703  setEquidistantU
704  (
705  U,
706  lenAcc,
707  maxIter,
708  spacingCorrInterval,
709  tolerance
710  );
711 
712  return U;
713 }
714 
715 
716 // * * * * * * * * * * * * * * Location Functions * * * * * * * * * * * * * //
717 
719 (
720  const scalar u,
721  const label CPI,
722  const label degree
723 ) const
724 {
725  return basis_.checkRange(u, CPI, degree);
726 }
727 
728 
729 bool NURBS3DCurve::checkRange(const scalar u, const label CPI) const
730 {
731  const label degree(basis_.degree());
732  return basis_.checkRange(u, CPI, degree);
733 }
734 
735 
736 scalar NURBS3DCurve::length(const label uIStart, const label uIEnd) const
737 {
738  // Compute derivatives wrt u for the given u interval.
739  const label lenSize(uIEnd - uIStart + 1);
740  vectorField dxdu(lenSize, Zero);
741  scalar length(Zero);
742 
743  forAll(dxdu, uI)
744  {
745  dxdu[uI] = curveDerivativeU(u_[uIStart + uI]);
746  }
747 
748  // Integrate.
749  for (label uI = 0; uI < (lenSize - 1); uI++)
750  {
751  length +=
752  0.5
753  *(mag(dxdu[uI + 1]) + mag(dxdu[uI]))
754  *(u_[uIStart + uI + 1]-u_[uIStart + uI]);
755  }
756 
757  return length;
758 }
759 
760 
762 (
763  const scalar uStart,
764  const scalar uEnd,
765  const label nPts
766 ) const
767 {
768  // Compute derivatives wrt u for the given u interval.
769  scalar length(Zero);
770  scalarField localU(nPts, Zero);
771  vectorField dxdu(nPts, Zero);
772 
773  forAll(localU, uI)
774  {
775  localU[uI] = uStart + scalar(uI)/scalar(nPts - 1)*(uEnd - uStart);
776  dxdu[uI] = curveDerivativeU(localU[uI]);
777  }
778 
779  // Integrate.
780  for (label uI = 0; uI < (nPts - 1); uI++)
781  {
782  length +=
783  0.5
784  *(mag(dxdu[uI + 1]) + mag(dxdu[uI]))
785  *(localU[uI + 1]-localU[uI]);
786  }
787 
788  return length;
789 }
790 
791 
792 scalar NURBS3DCurve::length() const
793 {
794  return length(scalar(0), (u_.size() - 1));
795 }
796 
797 
798 // * * * * * * * * * * * * * Derivative Functions * * * * * * * * * * * * * //
799 
801 {
802  const label degree(basis_.degree());
803  vector NWP(Zero);
804  vector dNduWP(Zero);
805  scalar NW(Zero);
806  scalar dNduW(Zero);
807 
808  forAll(CPs_, CPI)
809  {
810  const scalar basisValue(basis_.basisValue(CPI, degree, u));
811  const scalar basisDeriv(basis_.basisDerivativeU(CPI, degree, u));
812 
813  NWP += basisValue * weights_[CPI] * CPs_[CPI];
814  dNduWP += basisDeriv * weights_[CPI] * CPs_[CPI];
815  NW += basisValue * weights_[CPI];
816  dNduW += basisDeriv * weights_[CPI];
817  }
818 
819  const vector uDerivative((dNduWP - NWP*dNduW/NW)/NW);
820 
821  return uDerivative;
822 }
823 
824 
826 {
827  const label degree(basis_.degree());
828  vector NWP(Zero);
829  vector dNduWP(Zero);
830  vector d2Ndu2WP(Zero);
831  scalar NW(Zero);
832  scalar dNduW(Zero);
833  scalar d2Ndu2W(Zero);
834 
835  forAll(CPs_, CPI)
836  {
837  const scalar basisValue(basis_.basisValue(CPI, degree, u));
838  const scalar basisDeriv(basis_.basisDerivativeU(CPI, degree, u));
839  const scalar basis2Deriv(basis_.basisDerivativeUU(CPI, degree, u));
840 
841  NWP += basisValue * weights_[CPI] * CPs_[CPI];
842  dNduWP += basisDeriv * weights_[CPI] * CPs_[CPI];
843  d2Ndu2WP += basis2Deriv * weights_[CPI] * CPs_[CPI];
844  NW += basisValue * weights_[CPI];
845  dNduW += basisDeriv * weights_[CPI];
846  d2Ndu2W += basis2Deriv * weights_[CPI];
847  }
848 
849  const vector uuDerivative
850  (
851  (
852  d2Ndu2WP
853  - scalar(2)*dNduWP*dNduW/NW
854  - NWP*d2Ndu2W/NW
855  + scalar(2)*NWP*dNduW*dNduW/NW/NW
856  ) / NW
857  );
858 
859  return uuDerivative;
860 }
861 
862 
864 (
865  const scalar u,
866  const label CPI
867 )
868 {
869  // compute denominator.
870  const label degree(basis_.degree());
871  scalar NW(Zero);
872 
873  forAll(CPs_, CPJ)
874  {
875  NW += basis_.basisValue(CPJ, degree, u) * weights_[CPJ];
876  }
877 
878  const scalar basisValueI(basis_.basisValue(CPI, degree, u));
879  const scalar CPDerivative(basisValueI * weights_[CPI] / NW);
880 
881  return CPDerivative;
882 }
883 
884 
886 (
887  const scalar u,
888  const label CPI
889 )
890 {
891  // Compute nominator, denominator.
892  const label degree(basis_.degree());
893  vector NWP(Zero);
894  scalar NW(Zero);
895 
896  forAll(CPs_, CPJ)
897  {
898  const scalar basisValue(basis_.basisValue(CPJ, degree, u));
899  NWP += basisValue * weights_[CPJ] * CPs_[CPJ];
900  NW += basisValue * weights_[CPJ];
901  }
902 
903  // Compute derivative.
904  const scalar basisValueI(basis_.basisValue(CPI, degree, u));
905  const vector WDerivative(basisValueI/NW * (CPs_[CPI] - NWP/NW));
906 
907  return WDerivative;
908 }
909 
910 
912 (
913  const scalar uStart,
914  const scalar uEnd,
915  const label nPts
916 )
917 {
918  // compute derivatives wrt u for the given u interval
919  vectorField dxdu(nPts, Zero);
920  vectorField d2xdu2(nPts, Zero);
921  scalarField localU(nPts, Zero);
922  scalar lDerivative(Zero);
923 
924  forAll(localU, uI)
925  {
926  scalar& uLocal(localU[uI]);
927  uLocal = uStart + scalar(uI)/scalar(nPts - 1)*(uEnd - uStart);
928  dxdu[uI] = curveDerivativeU(uLocal);
929  d2xdu2[uI] = curveDerivativeUU(uLocal);
930  }
931 
932  // Integrate.
933  for (label uI = 0; uI<(nPts - 1); uI++)
934  {
935  lDerivative +=
936  0.5
937  * (
938  (dxdu[uI + 1]&d2xdu2[uI + 1])/mag(dxdu[uI + 1])
939  + (dxdu[uI]&d2xdu2[uI])/mag(dxdu[uI])
940  )
941  * (localU[uI + 1]-localU[uI]);
942  }
943 
944  return lDerivative;
945 }
946 
947 
948 // * * * * * * * * * * * * * * * Access Functions * * * * * * * * * * * * * //
949 
950 // -- Inlined in H--
951 
952 // * * * * * * * * * * * * * * * Write Functions * * * * * * * * * * * * * //
953 
955 {
956  write(name_);
957 }
958 
959 
961 (
962  const word fileName
963 )
964 {
965  if (Pstream::master())
966  {
967  OFstream curveFile(fileName);
968  OFstream curveFileCPs(fileName + "CPs");
969  OFstream curveFileBases(fileName + "Bases");
970  label degree(basis_.degree());
971 
972  vectorField& field(*this);
973 
974  forAll(*this, uI)
975  {
976  curveFile << field[uI].x() << " "
977  << field[uI].y() << " "
978  << field[uI].z()
979  << endl;
980  }
981 
982  forAll(CPs_, CPI)
983  {
984  curveFileCPs << CPs_[CPI].x() << " "
985  << CPs_[CPI].y() << " "
986  << CPs_[CPI].z()
987  << endl;
988  }
989 
990  forAll(*this, uI)
991  {
992  const scalar u(u_[uI]);
993  scalarField basesValues(CPs_.size());
994 
995  curveFileBases << u << " ";
996 
997  forAll(CPs_, CPI)
998  {
999  basesValues[CPI] = basis_.basisValue(CPI, degree, u);
1000  curveFileBases << basesValues[CPI] << " ";
1001  }
1002 
1003  curveFileBases << endl;
1004  }
1005  }
1006 }
1007 
1008 
1011  const fileName dirName,
1012  const fileName fileName
1013 )
1014 {
1015  if (Pstream::master())
1016  {
1017  OFstream curveFile(dirName/fileName);
1018  OFstream curveFileCPs(dirName/fileName + "CPs");
1019  OFstream curveFileBases(dirName/fileName + "Bases");
1020  label degree(basis_.degree());
1021 
1022  vectorField& field(*this);
1023 
1024  forAll(*this, uI)
1025  {
1026  curveFile << field[uI].x() << " "
1027  << field[uI].y() << " "
1028  << field[uI].z()
1029  << endl;
1030  }
1031 
1032  forAll(CPs_, CPI)
1033  {
1034  curveFileCPs << CPs_[CPI].x() << " "
1035  << CPs_[CPI].y() << " "
1036  << CPs_[CPI].z()
1037  << endl;
1038  }
1039 
1040  forAll(*this, uI)
1041  {
1042  const scalar u(u_[uI]);
1043  scalarField basesValues(CPs_.size());
1044 
1045  curveFileBases << u << " ";
1046 
1047  forAll(CPs_, CPI)
1048  {
1049  basesValues[CPI] = basis_.basisValue(CPI, degree, u);
1050  curveFileBases << basesValues[CPI] << " ";
1051  }
1052 
1053  curveFileBases << endl;
1054  }
1055  }
1056 }
1057 
1058 
1060 {
1061  writeWParses(name_);
1062 }
1063 
1064 
1067  const word fileName
1068 )
1069 {
1070  if (Pstream::master())
1071  {
1072  OFstream curveFile(fileName);
1073  OFstream curveFileCPs(fileName + "CPs");
1074  OFstream curveFileBases(fileName + "Bases");
1075  label degree(basis_.degree());
1076 
1077  vectorField& field(*this);
1078 
1079  forAll(*this, uI)
1080  {
1081  curveFile << "("
1082  << field[uI].x() << " "
1083  << field[uI].y() << " "
1084  << field[uI].z() << ")"
1085  << endl;
1086  }
1087 
1088  forAll(CPs_, CPI)
1089  {
1090  curveFileCPs << "("
1091  << CPs_[CPI].x() << " "
1092  << CPs_[CPI].y() << " "
1093  << CPs_[CPI].z() << ")"
1094  << endl;
1095  }
1096 
1097  forAll(*this, uI)
1098  {
1099  const scalar u(u_[uI]);
1100  scalarField basesValues(CPs_.size());
1101 
1102  curveFileBases << u << " ";
1103 
1104  forAll(CPs_, CPI)
1105  {
1106  basesValues[CPI] = basis_.basisValue(CPI, degree, u);
1107  curveFileBases << basesValues[CPI] << " ";
1108  }
1109 
1110  curveFileBases << endl;
1111  }
1112  }
1113 }
1114 
1115 
1118  const fileName dirName,
1119  const fileName fileName
1120 )
1121 {
1122  if (Pstream::master())
1123  {
1124  OFstream curveFile(dirName/fileName);
1125  OFstream curveFileCPs(dirName/fileName + "CPs");
1126  OFstream curveFileBases(dirName/fileName + "Bases");
1127  label degree(basis_.degree());
1128 
1129  vectorField& field(*this);
1130 
1131  forAll(*this, uI)
1132  {
1133  curveFile << "("
1134  << field[uI].x() << " "
1135  << field[uI].y() << " "
1136  << field[uI].z() << ")"
1137  << endl;
1138  }
1139 
1140  forAll(CPs_, CPI)
1141  {
1142  curveFileCPs << "("
1143  << CPs_[CPI].x() << " "
1144  << CPs_[CPI].y() << " "
1145  << CPs_[CPI].z() << ")"
1146  << endl;
1147  }
1148 
1149  forAll(*this, uI)
1150  {
1151  const scalar u(u_[uI]);
1152  scalarField basesValues(CPs_.size());
1153 
1154  curveFileBases << u << " ";
1155 
1156  forAll(CPs_, CPI)
1157  {
1158  basesValues[CPI] = basis_.basisValue(CPI, degree, u);
1159  curveFileBases << basesValues[CPI] << " ";
1160  }
1161 
1162  curveFileBases << endl;
1163  }
1164  }
1165 }
1166 
1167 
1168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1169 
1170 } // End namespace Foam
1171 
1172 // ************************************************************************* //
Foam::NURBS3DCurve::buildCurve
void buildCurve()
Build curve.
Definition: NURBS3DCurve.C:377
Foam::NURBS3DCurve::lengthDerivativeU
scalar lengthDerivativeU(const scalar uStart, const scalar uEnd, const label nPts)
Definition: NURBS3DCurve.C:912
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Foam::NURBSbasis::basisDerivativeU
scalar basisDerivativeU(const label iCP, const label degree, const scalar u) const
Basis derivative w.r.t u.
Definition: NURBSbasis.C:190
Foam::NURBS3DCurve::curveDerivativeUU
vector curveDerivativeUU(const scalar u) const
Curve second derivative wrt u at point ui.
Definition: NURBS3DCurve.C:825
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::NURBS3DCurve::nrm2D
const vector nrm2D(const scalar zVal, const scalar u) const
Definition: NURBS3DCurve.C:631
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::NURBS3DCurve::genEquidistant
scalarList genEquidistant(const label nPts=100, const label lenAcc=25, const label maxIter=10, const label spacingCorrInterval=-1, const scalar tolerance=1.e-5)
Definition: NURBS3DCurve.C:688
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::NURBS3DCurve::NURBS3DCurve
NURBS3DCurve(const NURBSbasis &basis, const List< vector > &CPs, const List< scalar > &weights, const scalarField &u, const label nPts, const word name="NURBS3DCurve")
Definition: NURBS3DCurve.C:210
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::NURBSbasis::basisValue
scalar basisValue(const label iCP, const label degree, const scalar u) const
Basis value.
Definition: NURBSbasis.C:140
Foam::NURBSbasis
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
Definition: NURBSbasis.H:55
Foam::NURBS3DCurve::setCPs
void setCPs(const List< vector > &CPs)
Set CPs.
Definition: NURBS3DCurve.C:359
Foam::NURBS3DCurve::nrm3D
const vector nrm3D(const vector &refTan, const scalar u) const
Find the normal to the curve, with the option of forcing a z-plane.
Definition: NURBS3DCurve.C:612
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::NURBSbasis::degree
const label & degree() const
Definition: NURBSbasisI.H:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::NURBS3DCurve::curveDerivativeU
vector curveDerivativeU(const scalar u) const
Curve derivative wrt u at point ui.
Definition: NURBS3DCurve.C:800
Foam::NURBS3DCurve::curvePoint
vector curvePoint(const scalar u) const
Curve point cartesian coordinates at ui.
Definition: NURBS3DCurve.C:457
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
NURBS3DCurve.H
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::NURBS3DCurve::setName
void setName(const word &name)
Set name.
Definition: NURBS3DCurve.C:371
Foam::NURBSbasis::basisDerivativeUU
scalar basisDerivativeUU(const label iCP, const label degree, const scalar u) const
Basis second derivative w.r.t u.
Definition: NURBSbasis.C:231
Foam::NURBS3DCurve::makeEquidistant
void makeEquidistant(const label lenAcc=25, const label maxIter=10, const label spacingCorrInterval=-1, const scalar tolerance=1.e-5)
Make curve points equidistant in cartesian space.
Definition: NURBS3DCurve.C:430
Foam::Field< scalar >
Foam::NURBS3DCurve::curveDerivativeCP
scalar curveDerivativeCP(const scalar u, const label CPI)
Definition: NURBS3DCurve.C:864
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::NURBS3DCurve::setNrm3DOrientation
void setNrm3DOrientation(const vector &givenNrm, const vector &givenTan)
Definition: NURBS3DCurve.C:291
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
field
rDeltaTY field()
Foam::NURBS3DCurve::write
void write()
Write curve to file.
Definition: NURBS3DCurve.C:954
Foam::NURBS3DCurve::findClosestCurvePoint
scalar findClosestCurvePoint(const vector &targetPoint, const label maxIter=1000, const scalar tolerance=1.e-13)
Definition: NURBS3DCurve.C:484
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::NURBS3DCurve::checkRange
bool checkRange(const scalar u, const label CPI, const label degree) const
Check if given parametric coordinate u and CP are linked.
Definition: NURBS3DCurve.C:719
Foam::NURBS3DCurve::setNrm2DOrientation
void setNrm2DOrientation(const vector &givenNrm, const scalar zVal)
Definition: NURBS3DCurve.C:317
U
U
Definition: pEqn.H:72
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
Foam::Vector< scalar >
Foam::NURBS3DCurve::length
scalar length() const
Calculate length for the entire curve length.
Definition: NURBS3DCurve.C:792
Foam::List< vector >
vectorField.H
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::curve
A single curve in a graph.
Definition: curve.H:58
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::NURBS3DCurve::setWeights
void setWeights(const List< scalar > &weights)
Set weights.
Definition: NURBS3DCurve.C:365
Foam::NURBS3DCurve::invert
void invert()
Invert CPs order and re-build curve.
Definition: NURBS3DCurve.C:408
Foam::NURBS3DCurve::writeWParses
void writeWParses()
Definition: NURBS3DCurve.C:1059
Foam::point
vector point
Point is a vector.
Definition: point.H:43
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::NURBSbasis::checkRange
bool checkRange(const scalar u, const label CPI, const label degree) const
Checks to see if given u is affected by given CP.
Definition: NURBSbasis.C:273
Foam::NURBS3DCurve::flipNrmOrientation
void flipNrmOrientation()
Flip the orientation of the nrm.
Definition: NURBS3DCurve.C:346
Foam::NURBS3DCurve::curveDerivativeWeight
vector curveDerivativeWeight(const scalar u, const label CPI)
Curve derivative wrt CPII at point u.
Definition: NURBS3DCurve.C:886
Foam::NURBS3DCurve::insertKnot
void insertKnot(const scalarField &oldKnots, const scalar uBar, const label kInsert)
Insert a knot by re-computing the control points.
Definition: NURBS3DCurve.C:646