Loading...
Searching...
No Matches
NURBS3DSurface.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-2022 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
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 "NURBS3DSurface.H"
31#include "vtkSurfaceWriter.H"
32#include "OFstream.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
41label NURBS3DSurface::sgn(const scalar val) const
42{
43 return val >= scalar(0) ? 1 : -1;
44}
45
46
47scalar NURBS3DSurface::abs(const scalar val) const
48{
49 return (sgn(val) == 1)? val: -val;
50}
51
52
53label NURBS3DSurface::mod(const label x, const label interval) const
54{
55 label ratio(x%interval);
56 return ratio < 0 ? ratio+interval : ratio;
57}
58
59
60void NURBS3DSurface::setCPUVLinking()
61{
62 const label uNCPs(uBasis_.nCPs());
63 const label vNCPs(vBasis_.nCPs());
64
65 CPsUCPIs_.setSize(uNCPs*vNCPs, -1);
66 CPsVCPIs_.setSize(uNCPs*vNCPs, -1);
67
68 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
69 {
70 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
71 {
72 const label CPI(vCPI*uNCPs + uCPI);
73 CPsUCPIs_[CPI] = uCPI;
74 CPsVCPIs_[CPI] = vCPI;
75 }
76 }
77}
78
79
80void NURBS3DSurface::setUniformUV
81(
82 scalarList& u,
83 scalarList& v,
84 const label nUPts,
85 const label nVPts
86) const
87{
88 u.setSize(nUPts*nVPts, Zero);
89 v.setSize(nUPts*nVPts, Zero);
90
91 for (label uI = 0; uI<nUPts; uI++)
92 {
93 scalar uValue = scalar(uI)/scalar(nUPts - 1);
94 for (label vI = 0; vI<nVPts; vI++)
95 {
96 const label ptI(uI*nVPts + vI);
97 u[ptI] = uValue;
98 v[ptI] = scalar(vI)/scalar(nVPts - 1);
99 }
100 }
101}
102
103
104void NURBS3DSurface::setUniformUV()
105{
106 setUniformUV(u_, v_, nUPts_, nVPts_);
107}
108
109
110bool NURBS3DSurface::bound
111(
112 scalar& u,
113 scalar& v,
114 const scalar minVal,
115 const scalar maxVal
116) const
117{
118 bool boundPoint =
119 boundDirection(u, minVal, maxVal)
120 || boundDirection(v, minVal, maxVal);
121
122 return boundPoint;
123}
124
125
126bool NURBS3DSurface::boundDirection
127(
128 scalar& u,
129 const scalar minVal,
130 const scalar maxVal
131) const
132{
133 bool boundPoint(false);
134
135 if (u < scalar(0))
136 {
137 u = minVal;
138 boundPoint = true;
139 }
140 else if (u > scalar(1))
141 {
142 u = maxVal;
143 boundPoint = true;
144 }
145
146 return boundPoint;
147}
148
149
150void NURBS3DSurface::setEquidistantR
151(
152 scalarList& R,
153 const scalar SHeld,
154 const label paramR,
155 const label lenAcc = 25,
156 const label maxIter = 10,
157 const label spacingCorrInterval = -1,
158 const scalar tolerance = 1.e-5
159) const
160{
161 const label nPts(R.size());
162 scalar SNull(SHeld);
163 scalar xLength(Zero);
164 const scalar rLength(scalar(1) / scalar(nPts - 1));
165
166 if (paramR == PARAMU)
167 {
168 xLength = lengthU(SHeld) / (nPts - 1);
169 }
170 else
171 {
172 xLength = lengthV(SHeld) / (nPts - 1);
173 }
174
175 R[0] = Zero;
176 R[nPts-1] = scalar(1);
177
178 for (label ptI=1; ptI<(nPts - 1); ptI++)
179 {
180 const scalar& RPrev(R[ptI - 1]);
181 scalar& RCurr(R[ptI]);
182 scalar direc(1);
183 scalar xDiff(0);
184 scalar delta(0);
185 bool overReached(false);
186
187 RCurr = RPrev + rLength;
188
189 // Find the starting U value to ensure target is within 1 rLength.
190 while (true)
191 {
192 bool bounded(false);
193
194 if (paramR == PARAMU)
195 {
196 bounded = bound(RCurr, SNull);
197 delta = lengthU(SHeld, RPrev, RCurr, lenAcc);
198 }
199 else
200 {
201 bounded = bound(SNull, RCurr);
202 delta = lengthV(SHeld, RPrev, RCurr, lenAcc);
203 }
204
205 xDiff = xLength - delta;
206
207 // Found the point.
208 if (abs(xDiff) < tolerance)
209 {
210 break;
211 }
212 else
213 {
214 direc = sgn(xDiff);
215
216 if (bounded && (direc == 1))
217 {
218 // rLength addition makes U exceed 1 so it becomes bounded.
219 // However, the desired x location still exceeds how far the
220 // bounded rLength can move (~e-5 error).
221 // Must force U to be u=0.999999.
222 overReached = true;
223 break;
224 }
225 else if (direc == scalar(1))
226 {
227 RCurr += rLength;
228 }
229 else
230 {
231 break;
232 }
233 }
234 }
235
236 if (!overReached)
237 {
238 label iter(0);
239
240 while (iter < maxIter)
241 {
242 // Set the new search length to ensure convergence and next v.
243 direc /= scalar(2);
244 RCurr += direc * rLength;
245
246 if (paramR == PARAMU)
247 {
248 bound(RCurr, SNull);
249 }
250 else
251 {
252 bound(SNull, RCurr);
253 }
254
255 // Can employ an occasional tolerance check from beg of curve.
256 if
257 (
258 (spacingCorrInterval != -1)
259 && (mod(ptI, spacingCorrInterval) == 0)
260 )
261 {
262 if (paramR == PARAMU)
263 {
264 delta = lengthU(SHeld, Zero, RCurr, ptI*lenAcc);
265 }
266 else
267 {
268 delta = lengthV(SHeld, Zero, RCurr, ptI*lenAcc);
269 }
270
271 xDiff = (ptI * xLength) - delta;
272 }
273 else
274 {
275 if (paramR == PARAMU)
276 {
277 delta = lengthU(SHeld, RPrev, RCurr, lenAcc);
278 }
279 else
280 {
281 delta = lengthV(SHeld, RPrev, RCurr, lenAcc);
282 }
283
284 xDiff = xLength - delta;
285 }
286
287 // Break if found point or find the new search direction.
288 if (abs(xDiff) < tolerance)
289 {
290 break;
291 }
292 else
293 {
294 const scalar oldDirec(direc);
295 direc = sgn(xDiff) * abs(oldDirec);
296 }
297
298 ++iter;
299 }
301 }
302}
303
304
305// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
306
308(
309 const List<vector>& CPs,
310 const label nUPts,
311 const label nVPts,
312 const NURBSbasis& uBasis,
313 const NURBSbasis& vBasis,
314 const word name
315)
316:
317 vectorField (nUPts*nVPts, Zero),
318
319 CPs_(CPs),
320 u_(nUPts*nVPts, Zero),
321 v_(nUPts*nVPts, Zero),
322 weights_(CPs.size(), scalar(1)),
323 nUPts_(nUPts),
324 nVPts_(nVPts),
325 name_(name),
326 uBasis_(uBasis),
327 vBasis_(vBasis),
328
329 givenInitNrm_(Zero),
330
331 CPsUCPIs_(0),
332 CPsVCPIs_(0),
333
334 nrmOrientation_(ALIGNED),
335
336 boundaryCPIDs_(nullptr)
338 setUniformUV();
339 buildSurface();
340 setCPUVLinking();
341}
342
343
345(
346 const List<vector>& CPs,
347 const List<scalar>& weights,
348 const label nUPts,
349 const label nVPts,
350 const NURBSbasis& uBasis,
351 const NURBSbasis& vBasis,
352 const word name
353)
354:
355 vectorField(nUPts*nVPts, Zero),
356
357 CPs_(CPs),
358 u_(nUPts*nVPts, Zero),
359 v_(nUPts*nVPts, Zero),
360 weights_(weights),
361 nUPts_(nUPts),
362 nVPts_(nVPts),
363 name_ (name),
364 uBasis_(uBasis),
365 vBasis_(vBasis),
366
367 givenInitNrm_(Zero),
368
369 CPsUCPIs_(0),
370 CPsVCPIs_(0),
371
372 nrmOrientation_(ALIGNED),
373
374 boundaryCPIDs_(nullptr)
376 setUniformUV();
377 buildSurface();
378 setCPUVLinking();
379}
380
381
383(
384 const List<vector>& CPs,
385 const label nUPts,
386 const label nVPts,
387 const label uDegree,
388 const label vDegree,
389 const label nCPsU,
390 const label nCPsV,
391 const word name
392)
393:
394 vectorField(nUPts*nVPts, Zero),
395
396 CPs_(CPs),
397 u_(nUPts*nVPts, Zero),
398 v_(nUPts*nVPts, Zero),
399 weights_(CPs.size(), scalar(1)),
400 nUPts_(nUPts),
401 nVPts_(nVPts),
402 name_(name),
403 uBasis_(nCPsU, uDegree),
404 vBasis_(nCPsV, vDegree),
405
406 givenInitNrm_(Zero),
407
408 CPsUCPIs_(0),
409 CPsVCPIs_(0),
410
411 nrmOrientation_(ALIGNED),
412
413 boundaryCPIDs_(nullptr)
414{
415 // Sanity checks
416 if (nCPsU*nCPsV != CPs_.size())
417 {
418 FatalErrorInFunction
419 << "nCPsU*nCPsV " << nCPsU*nCPsV
420 << " not equal to size of CPs " << CPs_.size()
421 << exit(FatalError);
422 }
423 // Initialize surface
424 setUniformUV();
425 buildSurface();
426 setCPUVLinking();
427}
428
429
431(
432 const List<vector>& CPs,
433 const label nUPts,
434 const label nVPts,
435 const label uDegree,
436 const label vDegree,
437 const label nCPsU,
438 const label nCPsV,
439 const scalarField& knotsU,
440 const scalarField& knotsV,
441 const word name
442)
443:
444 vectorField(nUPts*nVPts, Zero),
445
446 CPs_(CPs),
447 u_(nUPts*nVPts, Zero),
448 v_(nUPts*nVPts, Zero),
449 weights_(CPs.size(), scalar(1)),
450 nUPts_(nUPts),
451 nVPts_(nVPts),
452 name_(name),
453 uBasis_(nCPsU, uDegree, knotsU),
454 vBasis_(nCPsV, vDegree, knotsV),
455
456 givenInitNrm_(Zero),
457
458 CPsUCPIs_(0),
459 CPsVCPIs_(0),
460
461 nrmOrientation_(ALIGNED),
462
463 boundaryCPIDs_(nullptr)
464{
465 // Sanity checks
466 if (nCPsU*nCPsV != CPs_.size())
467 {
468 FatalErrorInFunction
469 << "nCPsU*nCPsV " << nCPsU*nCPsV
470 << " not equal to size of CPs " << CPs_.size()
471 << exit(FatalError);
472 }
473 // initialize surface
474 setUniformUV();
475 buildSurface();
476 setCPUVLinking();
477}
478
479
481(
482 const List<vector>& CPs,
483 const List<scalar>& weights,
484 const label nUPts,
485 const label nVPts,
486 const label uDegree,
487 const label vDegree,
488 const label nCPsU,
489 const label nCPsV,
490 const word name
491)
492:
493 vectorField(nUPts*nVPts, Zero),
494
495 CPs_(CPs),
496 u_(nUPts*nVPts, Zero),
497 v_(nUPts*nVPts, Zero),
498 weights_(weights),
499 nUPts_(nUPts),
500 nVPts_(nVPts),
501 name_(name),
502 uBasis_(nCPsU, uDegree),
503 vBasis_(nCPsV, vDegree),
504
505 givenInitNrm_(Zero),
506
507 CPsUCPIs_(0),
508 CPsVCPIs_(0),
509
510 nrmOrientation_(ALIGNED),
511
512 boundaryCPIDs_(nullptr)
513{
514 // Sanity checks
515 if (nCPsU*nCPsV != CPs_.size())
516 {
517 FatalErrorInFunction
518 << "nCPsU*nCPsV " << nCPsU*nCPsV
519 << " not equal to size of CPs " << CPs_.size()
520 << exit(FatalError);
521 }
522
523 // Initialize surface
524 setUniformUV();
525 buildSurface();
526 setCPUVLinking();
527}
528
529
531(
532 const List<vector>& CPs,
533 const List<scalar>& weights,
534 const label nUPts,
535 const label nVPts,
536 const label uDegree,
537 const label vDegree,
538 const label nCPsU,
539 const label nCPsV,
540 const scalarField& knotsU,
541 const scalarField& knotsV,
542 const word name
543)
544:
545 vectorField(nUPts*nVPts, Zero),
546
547 CPs_(CPs),
548 u_(nUPts*nVPts, Zero),
549 v_(nUPts*nVPts, Zero),
550 weights_(weights),
551 nUPts_(nUPts),
552 nVPts_(nVPts),
553 name_(name),
554 uBasis_(nCPsU, uDegree, knotsU),
555 vBasis_(nCPsV, vDegree, knotsV),
556
557 givenInitNrm_(Zero),
558
559 CPsUCPIs_(0),
560 CPsVCPIs_(0),
561
562 nrmOrientation_(ALIGNED),
563
564 boundaryCPIDs_(nullptr)
565{
566 // Sanity checks
567 if (nCPsU*nCPsV != CPs_.size())
568 {
569 FatalErrorInFunction
570 << "nCPsU*nCPsV " << nCPsU*nCPsV
571 << " not equal to size of CPs " << CPs_.size()
572 << exit(FatalError);
573 }
574 // Initialize surface
575 setUniformUV();
576 buildSurface();
577 setCPUVLinking();
578}
579
580// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
581
582// * * * * * * * * * * * * * * * * Set Functions * * * * * * * * * * * * * * //
583
585(
586 const vector& givenNrm,
587 const scalar u,
588 const scalar v
589)
590{
591 vector surfaceNrm(surfaceDerivativeU(u, v) ^ surfaceDerivativeV(u, v));
592
593 givenInitNrm_ = givenNrm;
594 surfaceNrm /= mag(surfaceNrm);
595
596 const scalar relation(givenNrm & surfaceNrm);
597
598 if (relation >= 0)
599 {
600 nrmOrientation_ = ALIGNED;
601 }
602 else
603 {
604 nrmOrientation_ = OPPOSED;
605 }
607 Info<< "Initial nrmOrientation after comparison to NURBS u="
608 << u << ",v=" << v << " nrm: " << nrmOrientation_
609 << endl;
610}
611
612
614{
615 if (nrmOrientation_ == ALIGNED)
616 {
617 nrmOrientation_ = OPPOSED;
618 }
619 else
620 {
621 nrmOrientation_ = ALIGNED;
622 }
623}
624
627{
628 CPs_ = CPs;
629}
630
632void NURBS3DSurface::setWeights(const scalarList& weights)
633{
634 weights_ = weights;
635}
636
639{
640 name_ = name;
641}
642
643
645{
646 const label uDegree(uBasis_.degree());
647 const label vDegree(vBasis_.degree());
648 const label uNCPs(uBasis_.nCPs());
649 const label vNCPs(vBasis_.nCPs());
650/*
651 Info<< "\nuDegree: " << uDegree << "\nvDegree: " << vDegree
652 << "\nnUPts: " << nUPts_ << "\nnVPts: " << nVPts_
653 << "\nuNCPs: " << uNCPs << "\nvNCPs: " << vNCPs
654 << "\nNURBSSurface:\nCPs: " << CPs
655 << endl;
656*/
657 vectorField& field = *this;
659
660 for (label uI = 0; uI<nUPts_; uI++)
661 {
662 for (label vI = 0; vI<nVPts_; vI++)
663 {
664 const label ptI(uI*nVPts_ + vI);
665 const scalar& u(u_[ptI]);
666 const scalar& v(v_[ptI]);
667 scalar NMW(Zero);
668
669 // Compute denominator.
670 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
671 {
672 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
673 {
674 const label CPI(vCPI*uNCPs + uCPI);
675
676 NMW +=
677 uBasis_.basisValue(uCPI, uDegree, u)
678 * vBasis_.basisValue(vCPI, vDegree, v)
679 * weights_[CPI];
680 }
681 }
682
683 // Compute the points.
684 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
685 {
686 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
687 {
688 const label CPI(vCPI*uNCPs + uCPI);
689
690 this->operator[](ptI) +=
691 CPs_[CPI]
692 * uBasis_.basisValue(uCPI, uDegree, u)
693 * vBasis_.basisValue(vCPI, vDegree, v)
694 * weights_[CPI]/NMW;
695 }
697 }
698 }
699}
700
701
702
704{
705 Info<< "Inverting NURBS surface " << name_ << " in u." << endl;
706
707 const label uNCPs(uBasis_.nCPs());
708 const label vNCPs(vBasis_.nCPs());
709 List<vector> invertedCPs(CPs_.size(), Zero);
710 List<scalar> invertedWeights(CPs_.size(), Zero);
711
712 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
713 {
714 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
715 {
716 const label CPI(vCPI*uNCPs + uCPI);
717 const label invUCPI(uNCPs-1-uCPI);
718 const label uInvCPI(vCPI*uNCPs + invUCPI);
719
720 invertedCPs[CPI] = CPs_[uInvCPI];
721 invertedWeights[CPI] = weights_[uInvCPI];
722 }
723 }
724
725 CPs_ = invertedCPs;
726 weights_ = invertedWeights;
727
728 buildSurface();
729}
730
731
733{
734 Info<< "Inverting NURBS surface " << name_ << " in v." << endl;
735
736 const label uNCPs(uBasis_.nCPs());
737 const label vNCPs(vBasis_.nCPs());
738 List<vector> invertedCPs(CPs_.size(), Zero);
739 List<scalar> invertedWeights(CPs_.size(), Zero);
740
741 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
742 {
743 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
744 {
745 const label CPI(vCPI*uNCPs + uCPI);
746 const label invVCPI(vNCPs-1-vCPI);
747 const label vInvCPI(invVCPI*uNCPs + uCPI);
748
749 invertedCPs[CPI] = CPs_[vInvCPI];
750 invertedWeights[CPI] = weights_[vInvCPI];
751 }
752 }
753
754 CPs_ = invertedCPs;
755 weights_ = invertedWeights;
756
757 buildSurface();
758}
759
760
762{
763 Info<< "Inverting NURBS surface " << name_ << " in u and v." << endl;
764
765 const label uNCPs(uBasis_.nCPs());
766 const label vNCPs(vBasis_.nCPs());
767 List<vector> invertedCPs(CPs_.size(), Zero);
768 List<scalar> invertedWeights(CPs_.size(), Zero);
769
770 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
771 {
772 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
773 {
774 const label CPI(vCPI*uNCPs + uCPI);
775 const label invUCPI(uNCPs - 1 - uCPI);
776 const label invVCPI(vNCPs - 1 - vCPI);
777 const label uvInvCPI(invVCPI*uNCPs + invUCPI);
778
779 invertedCPs[CPI] = CPs_[uvInvCPI];
780 invertedWeights[CPI] = weights_[uvInvCPI];
781 }
782 }
783
784 CPs_ = invertedCPs;
785 weights_ = invertedWeights;
786
787 buildSurface();
788}
789
790
792(
793 const label lenAcc,
794 const label maxIter,
795 const label spacingCorrInterval,
796 const scalar tolerance
797)
798{
799/*
800 Info<< "Making points equidistant is physical space on surface "
801 << name_
802 << endl;
803*/
804 // Equidistant spacing in u along v isoLines.
805 for (label vI = 0; vI<nVPts_; vI++)
806 {
807 scalarList UI(nUPts_, Zero);
808 const scalar VHeld(v_[vI]);
809 labelList uAddressing(nUPts_, -1);
810
811 // Set the point uAddressing to re-assign correct u_ values later.
812 forAll(uAddressing, uI)
813 {
814 const label ptI(uI*nVPts_ + vI);
815 uAddressing[uI] = ptI;
816 }
817 // Set equidistant u values.
818 setEquidistantR
819 (
820 UI,
821 VHeld,
822 PARAMU,
823 lenAcc,
824 maxIter,
825 spacingCorrInterval,
826 tolerance
827 );
828
829 // Re-assign new equidistant u values.
830 forAll(UI, uI)
831 {
832 const label& uAddress(uAddressing[uI]);
833 u_[uAddress] = UI[uI];
834 }
835 }
836
837 // Equidistant spacing in v along u isoLines.
838 for (label uI = 0; uI<nUPts_; uI++)
839 {
840 scalarList VI(nVPts_, Zero);
841 const scalar UHeld(u_[uI*nVPts_]);
842 labelList vAddressing(nUPts_, -1);
843
844 // Set the point vAddressing to re-assign correct u_ values later.
845 forAll(vAddressing, vI)
846 {
847 const label ptI(uI*nVPts_ + vI);
848 vAddressing[vI] = ptI;
849 }
850
851 // Set equidistant u values.
852 setEquidistantR
853 (
854 VI,
855 UHeld,
856 PARAMV,
857 lenAcc,
858 maxIter,
859 spacingCorrInterval,
860 tolerance
861 );
862
863 // Re-assign new equidistant u values.
864 forAll(VI, vI)
865 {
866 const label& vAddress(vAddressing[vI]);
867 v_[vAddress] = VI[vI];
868 }
869 }
871 buildSurface();
872}
873
874
875// * * * * * * * * * * * * * Point Calc Functions * * * * * * * * * * * * * //
876
878(
879 const scalar& u,
880 const scalar& v
881)
882{
883 const label uDegree(uBasis_.degree());
884 const label vDegree(vBasis_.degree());
885 const label uNCPs(uBasis_.nCPs());
886 const label vNCPs(vBasis_.nCPs());
887 scalar NMW(Zero);
888
889 // Compute denominator.
890 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
891 {
892 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
893 {
894 const label CPI(vCPI*uNCPs + uCPI);
895
896 NMW +=
897 uBasis_.basisValue(uCPI, uDegree, u)
898 * vBasis_.basisValue(vCPI, vDegree, v)
899 * weights_[CPI];
900 }
901 }
902
903 // Compute the points.
905
906 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
907 {
908 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
909 {
910 const label CPI(vCPI*uNCPs + uCPI);
911
912 point +=
913 CPs_[CPI]
914 * uBasis_.basisValue(uCPI, uDegree, u)
915 * vBasis_.basisValue(vCPI, vDegree, v)
916 * weights_[CPI]/NMW;
918 }
919
920 return point;
921}
922
923
925(
926 const vector& targetPoint,
927 const label maxIter,
928 const scalar tolerance
929)
930{
931 // Loop through surface points to find the closest one for initialization.
932 const vectorField& surface(*this);
933 scalar dist(GREAT);
934 label closePtI(-1);
935
936 forAll(surface, ptI)
937 {
938 const scalar distLoc(mag(targetPoint-surface[ptI]));
939
940 if (distLoc < dist)
941 {
942 dist = distLoc;
943 closePtI = ptI;
944 }
945 }
946
947 label iter(0);
948 scalar u(u_[closePtI]);
949 scalar v(v_[closePtI]);
950 vector xuv(surfacePoint(u, v));
951 scalar res(GREAT);
952 scalar resOld(GREAT);
953 scalar resDeriv(GREAT);
954 label nBoundsU(0);
955 label nBoundsV(0);
956
957 do
958 {
959 /*
960 const vector dxdu(surfaceDerivativeU(u, v));
961 const vector dxdv(surfaceDerivativeV(u, v));
962 const vector d2xdu2(surfaceDerivativeUU(u, v));
963 const vector d2xdv2(surfaceDerivativeVV(u, v));
964 const scalar uLHS((dxdu&dxdu) + ((xuv-targetPoint) & d2xdu2));
965 const scalar uRHS(-((xuv-targetPoint) & dxdu));
966 const scalar vLHS((dxdv&dxdv) + ((xuv-targetPoint) & d2xdv2));
967 const scalar vRHS(-((xuv-targetPoint) & dxdv));
968
969 // Update parametric coordinate and compute new point and
970 // bound param coordinates if needed.
971 // Compute residual.
972 u += uRHS/(uLHS+SMALL);
973 v += vRHS/(vLHS+SMALL);
974
975 bound(u, v);
976
977 xuv = surfacePoint(u, v);
978 res = mag((xuv-targetPoint) & surfaceDerivativeU(u, v))
979 + mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
980 */
981
982 const vector dxdu(surfaceDerivativeU(u, v));
983 const vector dxdv(surfaceDerivativeV(u, v));
984 const vector d2xdu2(surfaceDerivativeUU(u, v));
985 const vector d2xdv2(surfaceDerivativeVV(u, v));
986 const vector d2xduv(surfaceDerivativeUV(u, v));
987 const scalar a((dxdu&dxdu) + ((xuv-targetPoint) & d2xdu2));
988 const scalar b((dxdu&dxdv) + ((xuv-targetPoint) & d2xduv));
989 const scalar c=b;
990 const scalar d((dxdv&dxdv) + ((xuv-targetPoint) & d2xdv2));
991 const scalar invDenom = 1./(a*d-b*c);
992
993 const scalar uRHS(-((xuv-targetPoint) & dxdu));
994 const scalar vRHS(-((xuv-targetPoint) & dxdv));
995
996 // Update parametric coordinate and compute new point and
997 // bound param coordinates if needed.
998 // Compute residual.
999 u += ( d*uRHS-b*vRHS)*invDenom;
1000 v += (-c*uRHS+a*vRHS)*invDenom;
1001
1002 if (boundDirection(u))
1003 {
1004 nBoundsU++;
1005 }
1006 if (boundDirection(v))
1007 {
1008 nBoundsV++;
1009 }
1010
1011 xuv = surfacePoint(u, v);
1012 // If manual assignment in u is required, deal only with the v eqn
1013 if (nBoundsU >= 5)
1014 {
1015 res = mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
1016 resDeriv = mag(res-resOld)/resOld;
1017 resOld = res;
1018// Info<< targetPoint << " " << res << endl;
1019 }
1020 // If manual assignment in v is required, deal only with the u eqn
1021 else if (nBoundsV >= 5)
1022 {
1023 res = mag((xuv-targetPoint) & surfaceDerivativeU(u, v));
1024 resDeriv = mag(res-resOld)/resOld;
1025 resOld = res;
1026// Info<< targetPoint << " " << res << endl;
1027 }
1028 else if (nBoundsU <= 5 && nBoundsV <= 5)
1029 {
1030 res = mag((xuv-targetPoint) & surfaceDerivativeU(u, v))
1031 + mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
1032 resDeriv = mag(res-resOld)/resOld;
1033 resOld = res;
1034 }
1035 else
1036 {
1038 << "More than 5 bounds in both the u and v directions!"
1039 << "Something seems weird" << nBoundsU << " " << nBoundsV
1040 << endl;
1041 res = mag((xuv-targetPoint) & surfaceDerivativeU(u, v))
1042 + mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
1043 resDeriv = mag(res-resOld)/resOld;
1044 resOld = res;
1045 }
1046 }
1047 while ((iter++ < maxIter) && (res > tolerance));
1048
1049 scalarList closestParameters(2, u);
1050 closestParameters[1] = v;
1051 // warning if method did not reach threshold
1052 if (iter > maxIter)
1053 {
1055 << "Finding surface point closest to " << targetPoint
1056 << " for surface " << name_ << " failed \n"
1057 << " Number of bounding operations in u,v "
1058 << nBoundsU << " " << nBoundsV << endl
1059 << " Residual value and derivative " << res << " " << resDeriv
1060 << endl << endl;
1061 closestParameters = -1;
1062 }
1063
1064 return closestParameters;
1065}
1066
1067
1069(
1070 const vectorField& targetPoints,
1071 const label maxIter,
1072 const scalar tolerance
1073)
1074{
1075 auto tparamCoors = tmp<vector2DField>::New(targetPoints.size(), Zero);
1076 auto& paramCoors = tparamCoors.ref();
1077
1078 const vectorField& surface(*this);
1079 label nBoundedPoints(0);
1080 scalar maxResidual(0);
1081 scalar maxResidualDeriv(0);
1082 forAll(paramCoors, pI)
1083 {
1084 const vector& targetPoint(targetPoints[pI]);
1085
1086 // Loop through surface points to find the closest one for
1087 // initialization. The initialization could possibly be done with a
1088 // geodesical Laplace. Potentially faster?
1089 scalar dist(GREAT);
1090 label closePtI(-1);
1091
1092 forAll(surface, ptI)
1093 {
1094 const scalar distLoc(mag(targetPoint - surface[ptI]));
1095
1096 if (distLoc < dist)
1097 {
1098 dist = distLoc;
1099 closePtI = ptI;
1100 }
1101 }
1102
1103 label iter(0);
1104 scalar u(u_[closePtI]);
1105 scalar v(v_[closePtI]);
1106 vector xuv(surfacePoint(u, v));
1107 scalar res(GREAT);
1108 scalar resOld(GREAT);
1109 scalar resDeriv(GREAT);
1110 label nBoundsU(0);
1111 label nBoundsV(0);
1112
1113 do
1114 {
1115 const vector dxdu(surfaceDerivativeU(u, v));
1116 const vector dxdv(surfaceDerivativeV(u, v));
1117 const vector d2xdu2(surfaceDerivativeUU(u, v));
1118 const vector d2xdv2(surfaceDerivativeVV(u, v));
1119 const vector d2xduv(surfaceDerivativeUV(u, v));
1120 const scalar a((dxdu&dxdu) + ((xuv-targetPoint) & d2xdu2));
1121 const scalar b((dxdu&dxdv) + ((xuv-targetPoint) & d2xduv));
1122 const scalar c=b;
1123 const scalar d((dxdv&dxdv) + ((xuv-targetPoint) & d2xdv2));
1124 const scalar invDenom = 1./(a*d-b*c);
1125
1126 const scalar uRHS(-((xuv-targetPoint) & dxdu));
1127 const scalar vRHS(-((xuv-targetPoint) & dxdv));
1128
1129 // Update parametric coordinate and compute new point and
1130 // bound param coordinates if needed.
1131 // Compute residual.
1132 u += ( d*uRHS-b*vRHS)*invDenom;
1133 v += (-c*uRHS+a*vRHS)*invDenom;
1134
1135 if (boundDirection(u))
1136 {
1137 nBoundsU++;
1138 }
1139 if (boundDirection(v))
1140 {
1141 nBoundsV++;
1142 }
1143
1144 xuv = surfacePoint(u, v);
1145 // If manual assignment in u is required, deal only with the v eqn
1146 if (nBoundsU >= 5)
1147 {
1148 res = mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
1149 resDeriv = mag(res-resOld)/resOld;
1150 resOld = res;
1151 }
1152 // If manual assignment in b is required, deal only with the u eqn
1153 else if (nBoundsV >= 5)
1154 {
1155 res = mag((xuv-targetPoint) & surfaceDerivativeU(u, v));
1156 resDeriv = mag(res-resOld)/resOld;
1157 resOld = res;
1158 }
1159 else if (nBoundsU<=5 && nBoundsV<=5)
1160 {
1161 res = mag((xuv-targetPoint) & surfaceDerivativeU(u, v))
1162 + mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
1163 resDeriv = mag(res-resOld)/resOld;
1164 resOld = res;
1165 }
1166 else
1167 {
1169 << "More than 5 bounding operations in both the u and v directions!"
1170 << "u direction " << nBoundsU << endl
1171 << "v direction " << nBoundsV << endl
1172 << endl;
1173 res =
1174 mag((xuv-targetPoint) & surfaceDerivativeU(u, v))
1175 + mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
1176 resDeriv = mag(res-resOld)/resOld;
1177 resOld = res;
1178 }
1179 }
1180 while ((iter++ < maxIter) && (res > tolerance));
1181
1182 // warning if method did not reach threshold
1183 if (iter > maxIter)
1184 {
1185 nBoundedPoints++;
1186 maxResidual = max(maxResidual, res);
1187 maxResidualDeriv = max(maxResidualDeriv, resDeriv);
1188 }
1189
1190 paramCoors[pI].x() = u;
1191 paramCoors[pI].y() = v;
1192 }
1193
1194 Info<< "Points that couldn't reach the residual limit:: "
1195 << returnReduce(nBoundedPoints, sumOp<label>()) << nl
1196 << "Max residual after reaching iterations limit:: "
1197 << returnReduce(maxResidual, maxOp<scalar>()) << nl
1198 << "Max residual derivative after reaching iterations limit:: "
1199 << returnReduce(maxResidualDeriv, maxOp<scalar>()) << nl
1200 << endl;
1201
1202 return tparamCoors;
1203}
1204
1205
1207(
1208 const vector& targetPoint,
1209 const scalar& uInitGuess,
1210 const scalar& vInitGuess,
1211 const label maxIter,
1212 const scalar tolerance
1213)
1214{
1215 // Loop through surface points to find the closest one for initialization.
1216 label iter(0);
1217 scalar u(uInitGuess);
1218 scalar v(vInitGuess);
1219 vector xuv(surfacePoint(u, v));
1220 scalar res(GREAT);
1221
1222 do
1223 {
1224 const vector dxdu(surfaceDerivativeU(u, v));
1225 const vector dxdv(surfaceDerivativeV(u, v));
1226 const vector d2xdu2(surfaceDerivativeUU(u, v));
1227 const vector d2xdv2(surfaceDerivativeVV(u, v));
1228 const scalar uLHS((dxdu&dxdu) + ((xuv-targetPoint) & d2xdu2));
1229 const scalar uRHS(-((xuv-targetPoint) & dxdu));
1230 const scalar vLHS((dxdv&dxdv) + ((xuv-targetPoint) & d2xdv2));
1231 const scalar vRHS(-((xuv-targetPoint) & dxdv));
1232
1233 // Update parametric coordinate and compute new point and
1234 // bound param coordinates if needed.
1235 // Compute residual.
1236 u += uRHS/(uLHS+SMALL);
1237 v += vRHS/(vLHS+SMALL);
1238
1239 bound(u, v);
1240
1241 xuv = surfacePoint(u, v);
1242 res =
1243 mag((xuv-targetPoint) & surfaceDerivativeU(u, v))
1244 + mag((xuv-targetPoint) & surfaceDerivativeV(u, v));
1245 }
1246 while ((iter++ < maxIter) && (res > tolerance));
1247
1248 // warning if method did not reach threshold
1249 if (iter > maxIter)
1250 {
1252 << "Finding surface point closest to " << targetPoint << " failed."
1253 << endl;
1254 }
1255
1256 scalarList closestParameters(2, u);
1257 closestParameters[1] = v;
1258
1259 return closestParameters;
1260}
1261
1262
1263const vector NURBS3DSurface::nrm(scalar u, scalar v)
1264{
1265 vector surfaceNrm(Zero);
1266
1267 if (nrmOrientation_ == ALIGNED)
1268 {
1269 surfaceNrm = surfaceDerivativeU(u, v) ^ surfaceDerivativeV(u, v);
1270 }
1271 else
1272 {
1273 surfaceNrm = surfaceDerivativeV(u, v) ^ surfaceDerivativeU(u, v);
1274 }
1276 surfaceNrm /= mag(surfaceNrm);
1277
1278 return surfaceNrm;
1279}
1280
1281
1283(
1284 const label nUPts,
1285 const label nVPts,
1286 const label lenAcc,
1287 const label maxIter,
1288 const label spacingCorrInterval,
1289 const scalar tolerance
1290)
1291{
1292/*
1293 Info<< "Generating points equidistant in physical space on surface "
1294 << name_
1295 << endl;
1296*/
1297 // Preset U and V with uniform values.
1298 List<scalarList> UV(NPARAMS, scalarList(0));
1299
1300 scalarList& U(UV[PARAMU]);
1301 scalarList& V(UV[PARAMV]);
1302
1303 setUniformUV(U, V, nUPts, nVPts);
1304
1305 // Equidistant spacing in u along v isoLines.
1306 for (label vI = 0; vI<nVPts; vI++)
1307 {
1308 scalarList UI(nUPts, Zero);
1309 const scalar VHeld(V[vI]);
1310 labelList uAddressing(nUPts, -1);
1311
1312 // Set the point uAddressing to re-assign correct U values later.
1313 forAll(uAddressing, uI)
1314 {
1315 const label ptI(uI*nVPts + vI);
1316 uAddressing[uI] = ptI;
1317 }
1318 // Set equidistant u values.
1319 setEquidistantR
1320 (
1321 UI,
1322 VHeld,
1323 PARAMU,
1324 lenAcc,
1325 maxIter,
1326 spacingCorrInterval,
1327 tolerance
1328 );
1329
1330 // Re-assign new equidistant u values.
1331 forAll(UI, uI)
1332 {
1333 const label& uAddress(uAddressing[uI]);
1334 U[uAddress] = UI[uI];
1335 }
1336 }
1337
1338 // Equidistant spacing in v along u isoLines.
1339 for (label uI = 0; uI<nUPts; uI++)
1340 {
1341 scalarList VI(nVPts, Zero);
1342 const scalar UHeld(U[uI*nVPts]);
1343 labelList vAddressing(nUPts, -1);
1344
1345 // Set the point vAddressing to re-assign correct V values later.
1346 forAll(vAddressing, vI)
1347 {
1348 const label ptI(uI*nVPts + vI);
1349 vAddressing[vI] = ptI;
1350 }
1351
1352 // Set equidistant u values.
1353 setEquidistantR
1354 (
1355 VI,
1356 UHeld,
1357 PARAMV,
1358 lenAcc,
1359 maxIter,
1360 spacingCorrInterval,
1361 tolerance
1362 );
1363
1364 // Re-assign new equidistant u values.
1365 forAll(VI, vI)
1366 {
1367 const label& vAddress(vAddressing[vI]);
1368 V[vAddress] = VI[vI];
1369 }
1370 }
1372 return UV;
1373}
1374
1375
1376// * * * * * * * * * * * * * * Location Functions * * * * * * * * * * * * * //
1377
1379(
1380 const scalar u,
1381 const label CPI,
1382 const label uDegree
1383) const
1385 const label uCPI(CPsUCPIs_[CPI]);
1386
1387 return uBasis_.checkRange(u, uCPI, uDegree);
1388}
1389
1390
1392(
1393 const scalar u,
1394 const label CPI
1395) const
1397 const label uDegree(uBasis_.degree());
1398
1399 return checkRangeU(u, CPI, uDegree);
1400}
1401
1402
1404(
1405 const scalar v,
1406 const label CPI,
1407 const label vDegree
1408) const
1410 const label vCPI(CPsVCPIs_[CPI]);
1411
1412 return vBasis_.checkRange(v, vCPI, vDegree);
1413}
1414
1415
1417(
1418 const scalar v,
1419 const label CPI
1420) const
1422 const label vDegree(vBasis_.degree());
1423
1424 return checkRangeV(v, CPI, vDegree);
1425}
1426
1427
1429(
1430 const scalar v,
1431 const scalar u,
1432 const label CPI,
1433 const label uDegree,
1434 const label vDegree
1435) const
1436{
1437 if (checkRangeU(u, CPI, uDegree) && checkRangeV(v, CPI, vDegree))
1438 {
1439 return true;
1440 }
1441
1442 return false;
1443}
1444
1445
1447(
1448 const scalar v,
1449 const scalar u,
1450 const label CPI
1451) const
1452{
1453 const label uDegree(uBasis_.degree());
1454 const label vDegree(vBasis_.degree());
1455
1456 return checkRangeUV(u, v, CPI, uDegree, vDegree);
1457}
1458
1459
1461(
1462 const label vIConst,
1463 const label uIStart,
1464 const label uIEnd
1465) const
1466{
1467 // Compute derivatives wrt u for the given u interval.
1468 const label uLenSize(uIEnd - uIStart + 1);
1469 vectorField dxdu(uLenSize, Zero);
1470 scalar uLength(Zero);
1471
1472 forAll(dxdu, uI)
1473 {
1474 const label ptI((uIStart+uI)*nVPts_ + vIConst);
1475 const label& u(u_[ptI]);
1476 const label& v(v_[ptI]);
1477
1478 dxdu[uI] = surfaceDerivativeU(u, v);
1479 }
1480
1481 // Integrate.
1482 for (label uI = 0; uI<(uLenSize - 1); uI++)
1483 {
1484 const label ptI((uIStart+uI)*nVPts_ + vIConst);
1485
1486 uLength +=
1487 0.5*(mag(dxdu[uI + 1]) + mag(dxdu[uI]))*(u_[ptI + 1]-u_[ptI]);
1488 }
1489
1490 return uLength;
1491}
1492
1493
1495(
1496 const scalar vConst,
1497 const scalar uStart,
1498 const scalar uEnd,
1499 const label nPts
1500) const
1501{
1502 // Compute derivatives wrt u for the given u interval.
1503 vectorField dxdu(nPts, Zero);
1504 scalarField localU(nPts, Zero);
1505 scalar uLength(Zero);
1506
1507 forAll(localU, uI)
1508 {
1509 scalar& uLocal(localU[uI]);
1510 uLocal = uStart + scalar(uI)/scalar(nPts-1)*(uEnd-uStart);
1511 dxdu[uI] = surfaceDerivativeU(uLocal, vConst);
1512 }
1513
1514 // Integrate.
1515 for (label uI = 0; uI<(nPts - 1); uI++)
1516 {
1517 uLength +=
1518 0.5*(mag(dxdu[uI + 1]) + mag(dxdu[uI]))*(localU[uI + 1]-localU[uI]);
1519 }
1520
1521 return uLength;
1522}
1523
1525scalar NURBS3DSurface::lengthU(const label vIConst) const
1526{
1527 return lengthU(vIConst, 0, (nUPts_ - 1));
1528}
1529
1531scalar NURBS3DSurface::lengthU(const scalar vConst) const
1532{
1533 return lengthU(vConst, scalar(0), scalar(1), 100);
1534}
1535
1536
1538(
1539 const label uIConst,
1540 const label vIStart,
1541 const label vIEnd
1542) const
1543{
1544 // Compute derivatives wrt v for the given v interval.
1545 const label vLenSize(vIEnd - vIStart + 1);
1546 vectorField dxdv(vLenSize, Zero);
1547 scalar vLength(Zero);
1548
1549 forAll(dxdv, vI)
1550 {
1551 const label ptI((uIConst)*nVPts_ + (vIStart+vI));
1552 const label& u(u_[ptI]);
1553 const label& v(v_[ptI]);
1554
1555 dxdv[vI] = surfaceDerivativeV(u, v);
1556 }
1557
1558 // Integrate.
1559 for (label vI = 0; vI<(vLenSize - 1); vI++)
1560 {
1561 const label ptI((uIConst)*nVPts_ + (vIStart + vI));
1562
1563 vLength +=
1564 0.5*(mag(dxdv[vI + 1]) + mag(dxdv[vI]))*(v_[ptI + 1] - v_[ptI]);
1565 }
1566
1567 return vLength;
1568}
1569
1570
1572(
1573 const scalar uConst,
1574 const scalar vStart,
1575 const scalar vEnd,
1576 const label nPts
1577) const
1578{
1579 // Compute derivatives wrt v for the given v interval.
1580 vectorField dxdv(nPts, Zero);
1581 scalarField localV(nPts, Zero);
1582 scalar vLength(Zero);
1583
1584 forAll(localV, vI)
1585 {
1586 scalar& vLocal(localV[vI]);
1587 vLocal = vStart + scalar(vI)/scalar(nPts - 1)*(vEnd - vStart);
1588 dxdv[vI] = surfaceDerivativeV(uConst, vLocal);
1589 }
1590
1591 // Integrate.
1592 for (label vI = 0; vI<(nPts - 1); vI++)
1593 {
1594 vLength +=
1595 0.5*(mag(dxdv[vI + 1]) + mag(dxdv[vI]))*(localV[vI + 1]-localV[vI]);
1596 }
1597
1598 return vLength;
1599}
1600
1602scalar NURBS3DSurface::lengthV(const label uIConst) const
1603{
1604 return lengthV(uIConst, 0, (nVPts_ - 1));
1605}
1606
1607
1608scalar NURBS3DSurface::lengthV(const scalar uConst) const
1610 return lengthV(uConst, scalar(0), scalar(1), 100);
1611}
1612
1613
1614// * * * * * * * * * * * * * Derivative Functions * * * * * * * * * * * * * //
1615
1617(
1618 const scalar uIn,
1619 const scalar vIn
1620) const
1621{
1622 const label uDegree(uBasis_.degree());
1623 const label vDegree(vBasis_.degree());
1624 const label uNCPs(uBasis_.nCPs());
1625 const label vNCPs(vBasis_.nCPs());
1626 vector NMWP(Zero);
1627 vector dNduMWP(Zero);
1628 scalar NMW(Zero);
1629 scalar dNduMW(Zero);
1630
1631 scalar u = uIn;
1632 scalar v = vIn;
1633 bound(u, v);
1634
1635 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
1636 {
1637 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
1638 {
1639 const label CPI(vCPI*uNCPs + uCPI);
1640 const scalar uBasisValue(uBasis_.basisValue(uCPI, uDegree, u));
1641 const scalar vBasisValue(vBasis_.basisValue(vCPI, vDegree, v));
1642 const scalar uBasisDeriv
1643 (uBasis_.basisDerivativeU(uCPI, uDegree, u));
1644
1645 NMWP += uBasisValue * vBasisValue * weights_[CPI] * CPs_[CPI];
1646 dNduMWP += uBasisDeriv * vBasisValue * weights_[CPI] * CPs_[CPI];
1647 NMW += uBasisValue * vBasisValue * weights_[CPI];
1648 dNduMW += uBasisDeriv * vBasisValue * weights_[CPI];
1649 }
1650 }
1652 const vector uDerivative((dNduMWP - dNduMW*NMWP/(NMW+SMALL))/(NMW+SMALL));
1653
1654 return uDerivative;
1655}
1656
1657
1659(
1660 const scalar uIn,
1661 const scalar vIn
1662) const
1663{
1664 const label uDegree(uBasis_.degree());
1665 const label vDegree(vBasis_.degree());
1666 const label uNCPs(uBasis_.nCPs());
1667 const label vNCPs(vBasis_.nCPs());
1668 vector NMWP(Zero);
1669 vector dMdvNWP(Zero);
1670 scalar NMW(Zero);
1671 scalar dMdvNW(Zero);
1672
1673 scalar u = uIn;
1674 scalar v = vIn;
1675 bound(u, v);
1676
1677 for (label vCPI = 0; vCPI<vNCPs; vCPI++)
1678 {
1679 for (label uCPI = 0; uCPI<uNCPs; uCPI++)
1680 {
1681 const label CPI(vCPI*uNCPs + uCPI);
1682 const scalar uBasisValue(uBasis_.basisValue(uCPI, uDegree, u));
1683 const scalar vBasisValue(vBasis_.basisValue(vCPI, vDegree, v));
1684 const scalar vBasisDeriv
1685 (vBasis_.basisDerivativeU(vCPI, vDegree, v));
1686
1687 NMWP += uBasisValue * vBasisValue * weights_[CPI] * CPs_[CPI];
1688 dMdvNWP += vBasisDeriv * uBasisValue * weights_[CPI] * CPs_[CPI];
1689 NMW += uBasisValue * vBasisValue * weights_[CPI];
1690 dMdvNW += vBasisDeriv * uBasisValue * weights_[CPI];
1691 }
1692 }
1694 const vector vDerivative((dMdvNWP - dMdvNW*NMWP/(NMW+SMALL))/(NMW+SMALL));
1695
1696 return vDerivative;
1697}
1698
1699
1701(
1702 const scalar uIn,
1703 const scalar vIn
1704) const
1705{
1706 const label uDegree(uBasis_.degree());
1707 const label vDegree(vBasis_.degree());
1708 const label uNCPs(uBasis_.nCPs());
1709 const label vNCPs(vBasis_.nCPs());
1710 vector NMWP(Zero);
1711 vector dNduMWP(Zero);
1712 vector dMdvNWP(Zero);
1713 vector dNMduvWP(Zero);
1714 scalar NMW(Zero);
1715 scalar dNduMW(Zero);
1716 scalar dMdvNW(Zero);
1717 scalar dNMduvW(Zero);
1718
1719 scalar u = uIn;
1720 scalar v = vIn;
1721 bound(u, v);
1722
1723 for (label vCPI=0; vCPI<vNCPs; vCPI++)
1724 {
1725 for (label uCPI=0; uCPI<uNCPs; uCPI++)
1726 {
1727 const label CPI(vCPI*uNCPs + uCPI);
1728 const scalar uBasisValue(uBasis_.basisValue(uCPI, uDegree, u));
1729 const scalar vBasisValue(vBasis_.basisValue(vCPI, vDegree, v));
1730 const scalar uBasisDeriv
1731 (uBasis_.basisDerivativeU(uCPI, uDegree, u));
1732 const scalar vBasisDeriv
1733 (vBasis_.basisDerivativeU(vCPI, vDegree, v));
1734 //Info<< "vCPI=" << vCPI << ",uCPI=" << uCPI
1735 // << " N=" << uBasisValue << ",N'=" << uBasisDeriv
1736 // << " M=" << vBasisValue << ",M'=" << vBasisDeriv
1737 // << endl;
1738 NMWP += vBasisValue * uBasisValue * weights_[CPI] * CPs_[CPI];
1739 dNduMWP += uBasisDeriv * vBasisValue * weights_[CPI] * CPs_[CPI];
1740 dMdvNWP += vBasisDeriv * uBasisValue * weights_[CPI] * CPs_[CPI];
1741 dNMduvWP += uBasisDeriv * vBasisDeriv * weights_[CPI] * CPs_[CPI];
1742 NMW += vBasisValue * uBasisValue * weights_[CPI];
1743 dNduMW += uBasisDeriv * vBasisValue * weights_[CPI];
1744 dMdvNW += vBasisDeriv * uBasisValue * weights_[CPI];
1745 dNMduvW += uBasisDeriv * vBasisDeriv * weights_[CPI];
1746 }
1747 }
1748
1749 const vector uvDerivative
1750 (
1751 (
1752 dNMduvWP
1753 - (dNMduvW*NMWP + dMdvNW*dNduMWP + dNduMW*dMdvNWP)/(NMW+SMALL)
1754 + scalar(2)*dNduMW*dMdvNW*NMWP/(NMW+SMALL)/(NMW+SMALL)
1755 ) / (NMW+SMALL)
1756 );
1757
1758 return uvDerivative;
1759}
1760
1761
1763(
1764 const scalar uIn,
1765 const scalar vIn
1766) const
1767{
1768 const label uDegree(uBasis_.degree());
1769 const label vDegree(vBasis_.degree());
1770 const label uNCPs(uBasis_.nCPs());
1771 const label vNCPs(vBasis_.nCPs());
1772 vector NMWP(Zero);
1773 vector dNduMWP(Zero);
1774 vector d2Ndu2MWP(Zero);
1775 scalar NMW(Zero);
1776 scalar dNduMW(Zero);
1777 scalar d2Ndu2MW(Zero);
1778
1779 scalar u = uIn;
1780 scalar v = vIn;
1781 bound(u, v);
1782
1783 for (label vCPI=0; vCPI<vNCPs; vCPI++)
1784 {
1785 for (label uCPI=0; uCPI<uNCPs; uCPI++)
1786 {
1787 const label CPI(vCPI*uNCPs + uCPI);
1788 const scalar uBasisValue(uBasis_.basisValue(uCPI, uDegree, u));
1789 const scalar vBasisValue(vBasis_.basisValue(vCPI, vDegree, v));
1790 const scalar uBasisDeriv
1791 (uBasis_.basisDerivativeU(uCPI, uDegree, u));
1792 const scalar uBasis2Deriv
1793 (uBasis_.basisDerivativeUU(uCPI, uDegree, u));
1794
1795 NMWP += uBasisValue * vBasisValue * weights_[CPI] * CPs_[CPI];
1796 dNduMWP += uBasisDeriv * vBasisValue * weights_[CPI] * CPs_[CPI];
1797 d2Ndu2MWP += uBasis2Deriv * vBasisValue * weights_[CPI] * CPs_[CPI];
1798 NMW += uBasisValue * vBasisValue * weights_[CPI];
1799 dNduMW += uBasisDeriv * vBasisValue * weights_[CPI];
1800 d2Ndu2MW += uBasis2Deriv * vBasisValue * weights_[CPI];
1801 }
1802 }
1803
1804 const vector uuDerivative
1805 (
1806 (
1807 d2Ndu2MWP
1808 - scalar(2)*dNduMW*dNduMWP/(NMW+SMALL)
1809 - d2Ndu2MW*NMWP/(NMW+SMALL)
1810 + scalar(2)*dNduMW*dNduMW*NMWP/(NMW+SMALL)/(NMW+SMALL)
1811 ) / (NMW+SMALL)
1812 );
1813
1814 return uuDerivative;
1815}
1816
1817
1819(
1820 const scalar uIn,
1821 const scalar vIn
1822) const
1823{
1824 const label uDegree(uBasis_.degree());
1825 const label vDegree(vBasis_.degree());
1826 const label uNCPs(uBasis_.nCPs());
1827 const label vNCPs(vBasis_.nCPs());
1828 vector NMWP(Zero);
1829 vector dMdvNWP(Zero);
1830 vector d2Mdv2NWP(Zero);
1831 scalar NMW(Zero);
1832 scalar dMdvNW(Zero);
1833 scalar d2Mdv2NW(Zero);
1834
1835 scalar u = uIn;
1836 scalar v = vIn;
1837 bound(u, v);
1838
1839 for (label vCPI=0; vCPI<vNCPs; vCPI++)
1840 {
1841 for (label uCPI=0; uCPI<uNCPs; uCPI++)
1842 {
1843 const label CPI(vCPI*uNCPs + uCPI);
1844 const scalar uBasisValue(uBasis_.basisValue(uCPI, uDegree, u));
1845 const scalar vBasisValue(vBasis_.basisValue(vCPI, vDegree, v));
1846 const scalar vBasisDeriv
1847 (vBasis_.basisDerivativeU(vCPI, vDegree, v));
1848 const scalar vBasis2Deriv
1849 (vBasis_.basisDerivativeUU(vCPI, vDegree, v));
1850
1851 NMWP += vBasisValue * uBasisValue * weights_[CPI] * CPs_[CPI];
1852 dMdvNWP += vBasisDeriv * uBasisValue * weights_[CPI] * CPs_[CPI];
1853 d2Mdv2NWP += vBasis2Deriv * uBasisValue * weights_[CPI] * CPs_[CPI];
1854 NMW += vBasisValue * uBasisValue * weights_[CPI];
1855 dMdvNW += vBasisDeriv * uBasisValue * weights_[CPI];
1856 d2Mdv2NW += vBasis2Deriv * uBasisValue * weights_[CPI];
1857 }
1858 }
1859
1860 const vector vvDerivative
1861 (
1862 (
1863 d2Mdv2NWP
1864 - scalar(2)*dMdvNW*dMdvNWP/(NMW+SMALL)
1865 - d2Mdv2NW*NMWP/(NMW+SMALL)
1866 + scalar(2)*dMdvNW*dMdvNW*NMWP/(NMW+SMALL)/(NMW+SMALL)
1867 ) / (NMW+SMALL)
1868 );
1869
1870 return vvDerivative;
1871}
1872
1873
1875(
1876 const scalar u,
1877 const scalar v,
1878 const label CPI
1879) const
1880{
1881 //Info<< "u,v,cpI " << u << " " << v << " " << CPI << endl;
1882 // compute denominator.
1883 const label uDegree(uBasis_.degree());
1884 const label vDegree(vBasis_.degree());
1885 const label uNCPs(uBasis_.nCPs());
1886 const label vNCPs(vBasis_.nCPs());
1887 const label uCPI(CPsUCPIs_[CPI]);
1888 const label vCPI(CPsVCPIs_[CPI]);
1889 scalar NMW(Zero);
1890
1891 for (label vCPJ=0; vCPJ<vNCPs; vCPJ++)
1892 {
1893 for (label uCPJ=0; uCPJ<uNCPs; uCPJ++)
1894 {
1895 const label CPJ(vCPJ*uNCPs + uCPJ);
1896 const scalar uBasisValue(uBasis_.basisValue(uCPJ, uDegree, u));
1897 const scalar vBasisValue(vBasis_.basisValue(vCPJ, vDegree, v));
1898
1899 NMW += vBasisValue * uBasisValue * weights_[CPJ];
1900 }
1901 }
1902 //Info<< "denom " << NMW << endl;
1903
1904 const scalar CPDerivative
1905 (
1906 uBasis_.basisValue(uCPI, uDegree, u)
1907 * vBasis_.basisValue(vCPI, vDegree, v)
1908 * weights_[CPI]
1909 / (NMW+SMALL)
1910 );
1911
1912 return CPDerivative;
1913}
1914
1915
1917(
1918 const scalar u,
1919 const scalar v,
1920 const label CPI
1921) const
1922{
1923 // Compute nominator, denominator.
1924 const label uDegree(uBasis_.degree());
1925 const label vDegree(vBasis_.degree());
1926 const label uNCPs(uBasis_.nCPs());
1927 const label vNCPs(vBasis_.nCPs());
1928 const label uCPI(CPsUCPIs_[CPI]);
1929 const label vCPI(CPsVCPIs_[CPI]);
1930 vector NMWP(Zero);
1931 scalar NMW(Zero);
1932
1933 for (label vCPJ=0; vCPJ<vNCPs; vCPJ++)
1934 {
1935 for (label uCPJ=0; uCPJ<uNCPs; uCPJ++)
1936 {
1937 const label CPJ(vCPJ*uNCPs + uCPJ);
1938 const scalar uBasisValue(uBasis_.basisValue(uCPJ, uDegree, u));
1939 const scalar vBasisValue(vBasis_.basisValue(vCPJ, vDegree, v));
1940
1941 NMWP += uBasisValue * vBasisValue * weights_[CPJ] * CPs_[CPJ];
1942 NMW += uBasisValue * vBasisValue * weights_[CPJ];
1943 }
1944 }
1945
1946 // Compute derivative.
1947 const vector WDerivative
1948 (
1949 uBasis_.basisValue(uCPI, uDegree, u)
1950 * vBasis_.basisValue(vCPI, vDegree, v)
1951 * (CPs_[CPI] - (NMWP / (NMW+SMALL)))
1952 / (NMW+SMALL)
1953 );
1954
1955 return WDerivative;
1956}
1957
1958
1960(
1961 const scalar vConst,
1962 const scalar uStart,
1963 const scalar uEnd,
1964 const label nPts
1965) const
1966{
1967 // compute derivatives wrt u for the given u interval
1968 vectorField dxdu(nPts, Zero);
1969 vectorField d2xdu2(nPts, Zero);
1970 scalarField localU(nPts, Zero);
1971 scalar ulDerivative(Zero);
1972
1973 forAll(localU, uI)
1974 {
1975 scalar& uLocal(localU[uI]);
1976 uLocal = uStart + scalar(uI)/scalar(nPts-1)*(uEnd-uStart);
1977 dxdu[uI] = surfaceDerivativeU(uLocal, vConst);
1978 d2xdu2[uI] = surfaceDerivativeUU(uLocal, vConst);
1979 }
1980
1981 // Integrate.
1982 for (label uI=0; uI<(nPts-1); uI++)
1983 {
1984 ulDerivative +=
1985 0.5
1986 * (
1987 (dxdu[uI+1]&d2xdu2[uI+1])/(mag(dxdu[uI+1])+SMALL)
1988 + (dxdu[uI ]&d2xdu2[uI ])/(mag(dxdu[uI ])+SMALL)
1989 )
1990 * (localU[uI+1]-localU[uI]);
1991 }
1992
1993 return ulDerivative;
1994}
1995
1996
1998(
1999 const scalar uConst,
2000 const scalar vStart,
2001 const scalar vEnd,
2002 const label nPts
2003) const
2004{
2005 // Compute derivatives wrt v for the given v interval.
2006 vectorField dxdv(nPts, Zero);
2007 vectorField d2xdv2(nPts, Zero);
2008 scalarField localV(nPts, Zero);
2009 scalar vlDerivative(Zero);
2010
2011 forAll(localV, vI)
2012 {
2013 scalar& vLocal(localV[vI]);
2014 vLocal = vStart + scalar(vI)/scalar(nPts-1)*(vEnd-vStart);
2015 dxdv[vI] = surfaceDerivativeV(uConst, vLocal);
2016 d2xdv2[vI] = surfaceDerivativeVV(uConst, vLocal);
2017 }
2018
2019 // Integrate.
2020 for (label vI=0; vI<(nPts-1); vI++)
2021 {
2022 vlDerivative +=
2023 0.5
2024 * (
2025 (dxdv[vI+1]&d2xdv2[vI+1])/(mag(dxdv[vI+1])+SMALL)
2026 + (dxdv[vI ]&d2xdv2[vI ])/(mag(dxdv[vI ])+SMALL)
2027 )
2028 * (localV[vI+1]-localV[vI]);
2029 }
2031 return vlDerivative;
2032}
2033
2034
2035// * * * * * * * * * * * * * * * Access Functions * * * * * * * * * * * * * //
2036
2038{
2039 if (!boundaryCPIDs_)
2040 {
2041 const label uNCPs(uBasis_.nCPs());
2042 const label vNCPs(vBasis_.nCPs());
2043 const label nBoundCPs(2*uNCPs+2*vNCPs-4);
2044 boundaryCPIDs_.reset(new labelList(nBoundCPs, -1));
2045 whichBoundaryCPID_.reset(new labelList(uNCPs*vNCPs, -1));
2046
2047 // v-constant cps
2048 label bID(0);
2049 for (label vI=0; vI<vNCPs; vI+=vNCPs-1)
2050 {
2051 for (label uI=0; uI<uNCPs; uI++)
2052 {
2053 const label CPI(vI*uNCPs + uI);
2054 whichBoundaryCPID_()[CPI] = bID;
2055 boundaryCPIDs_()[bID++] = CPI;
2056 }
2057 }
2058 // u-constant cps
2059 for (label uI=0; uI<uNCPs; uI+=uNCPs-1)
2060 {
2061 // corner CPS already accounted for
2062 for (label vI=1; vI<vNCPs-1; vI++)
2063 {
2064 const label CPI(vI*uNCPs + uI);
2065 whichBoundaryCPID_()[CPI] = bID;
2066 boundaryCPIDs_()[bID++] = CPI;
2067 }
2069 }
2070
2071 return boundaryCPIDs_();
2072}
2073
2076{
2077 return getBoundaryCPIDs();
2078}
2079
2080
2081const label& NURBS3DSurface::whichBoundaryCPI(const label& globalCPI)
2082{
2083 if (!whichBoundaryCPID_)
2084 {
2086 }
2088 return whichBoundaryCPID_()[globalCPI];
2089}
2090
2091
2092// * * * * * * * * * * * * * * * Write Functions * * * * * * * * * * * * * //
2095{
2096 write(name_);
2097}
2098
2099
2101{
2102 if (Pstream::master())
2103 {
2104 OFstream surfaceFile(fileName);
2105 OFstream surfaceFileCPs(fileName+"CPs");
2106 vectorField& surface(*this);
2107
2108 forAll(*this, ptI)
2109 {
2110 surfaceFile
2111 << surface[ptI].x() << " "
2112 << surface[ptI].y() << " "
2113 << surface[ptI].z()
2114 << endl;
2115 }
2116
2117 forAll(CPs_, CPI)
2118 {
2119 surfaceFileCPs
2120 << CPs_[CPI].x() << " "
2121 << CPs_[CPI].y() << " "
2122 << CPs_[CPI].z()
2123 << endl;
2124 }
2125/*
2126 const label uDegree(uBasis_.degree());
2127 const label vDegree(vBasis_.degree());
2128 const label uNCPs(uBasis_.nCPs());
2129 const label vNCPs(vBasis_.nCPs());
2130
2131 OFstream surfaceFileUBases(fileName+"UBases");
2132 OFstream surfaceFileVBases(fileName+"VBases");
2133
2134 forAll(*this, ptI)
2135 {
2136 const scalar& u(u_[ptI]);
2137 const scalar& v(v_[ptI]);
2138 scalarField uBasesValues(uNCPs);
2139 scalarField vBasesValues(vNCPs);
2140
2141 surfaceFileUBases << u << " ";
2142 surfaceFileVBases << v << " ";
2143
2144 forAll(CPs_, CPI)
2145 {
2146 basesValues[CPI] = basis_.basisValue(CPI, degree, u);
2147 curveFileBases << basesValues[CPI] << " ";
2148 }
2149
2150 curveFileBases << endl;
2151 }
2152*/
2153 }
2154}
2155
2156
2157void NURBS3DSurface::write(const fileName dirName, const fileName fileName)
2158{
2159 if (Pstream::master())
2160 {
2161 OFstream surfaceFile(dirName/fileName);
2162 OFstream surfaceFileCPs(dirName/fileName+"CPs");
2163 vectorField& surface(*this);
2164
2165 forAll(*this, ptI)
2166 {
2167 surfaceFile
2168 << surface[ptI].x() << " "
2169 << surface[ptI].y() << " "
2170 << surface[ptI].z()
2171 << endl;
2172 }
2173
2174 forAll(CPs_, CPI)
2175 {
2176 surfaceFileCPs
2177 << CPs_[CPI].x() << " "
2178 << CPs_[CPI].y() << " "
2179 << CPs_[CPI].z()
2180 << endl;
2181 }
2182/*
2183 const label uDegree(uBasis_.degree());
2184 const label vDegree(vBasis_.degree());
2185 const label uNCPs(uBasis_.nCPs());
2186 const label vNCPs(vBasis_.nCPs());
2187
2188 OFstream surfaceFileUBases(dirName/fileName+"UBases");
2189 OFstream surfaceFileVBases(dirName/fileName+"VBases");
2190
2191 forAll(*this, ptI)
2192 {
2193 const scalar& u(u_[ptI]);
2194 const scalar& v(v_[ptI]);
2195 scalarField uBasesValues(uNCPs);
2196 scalarField vBasesValues(vNCPs);
2197
2198 surfaceFileUBases << u << " ";
2199 surfaceFileVBases << v << " ";
2200
2201 forAll(CPs_, CPI)
2202 {
2203 basesValues[CPI] = basis_.basisValue(CPI, degree, u);
2204 curveFileBases << basesValues[CPI] << " ";
2205 }
2206
2207 curveFileBases << endl;
2208 }
2209*/
2210 }
2211}
2212
2215{
2216 writeWParses(name_);
2217}
2218
2219
2221{
2222 if (Pstream::master())
2223 {
2224 OFstream surfaceFile(fileName);
2225 OFstream surfaceFileCPs(fileName+"CPs");
2226 vectorField& surface(*this);
2227
2228 forAll(*this, ptI)
2229 {
2230 surfaceFile
2231 << "("
2232 << surface[ptI].x() << " "
2233 << surface[ptI].y() << " "
2234 << surface[ptI].z() << ")"
2235 << endl;
2236 }
2237
2238 forAll(CPs_, CPI)
2239 {
2240 surfaceFileCPs
2241 << "("
2242 << CPs_[CPI].x() << " "
2243 << CPs_[CPI].y() << " "
2244 << CPs_[CPI].z() << ")"
2245 << endl;
2246 }
2247/*
2248 const label uDegree(uBasis_.degree());
2249 const label vDegree(vBasis_.degree());
2250 const label uNCPs(uBasis_.nCPs());
2251 const label vNCPs(vBasis_.nCPs());
2252
2253 OFstream surfaceFileUBases(fileName+"UBases");
2254 OFstream surfaceFileVBases(fileName+"VBases");
2255
2256 forAll(*this, ptI)
2257 {
2258 const scalar& u(u_[ptI]);
2259 const scalar& v(v_[ptI]);
2260 scalarField uBasesValues(uNCPs);
2261 scalarField vBasesValues(vNCPs);
2262
2263 surfaceFileUBases << u << " ";
2264 surfaceFileVBases << v << " ";
2265
2266 forAll(CPs_, CPI)
2267 {
2268 basesValues[CPI] = basis_.basisValue(CPI, degree, u);
2269 curveFileBases << basesValues[CPI] << " ";
2270 }
2271
2272 curveFileBases << endl;
2273 }
2274*/
2275 }
2276}
2277
2278
2280(
2281 const fileName dirName,
2282 const fileName fileName
2283)
2284{
2285 if (Pstream::master())
2286 {
2287 OFstream surfaceFile(dirName/fileName);
2288 OFstream surfaceFileCPs(dirName/fileName+"CPs");
2289 vectorField& surface(*this);
2290
2291 forAll(*this, ptI)
2292 {
2293 surfaceFile
2294 << "("
2295 << surface[ptI].x() << " "
2296 << surface[ptI].y() << " "
2297 << surface[ptI].z() << ")"
2298 << endl;
2299 }
2300
2301 forAll(CPs_, CPI)
2302 {
2303 surfaceFileCPs
2304 << "("
2305 << CPs_[CPI].x() << " "
2306 << CPs_[CPI].y() << " "
2307 << CPs_[CPI].z() << ")"
2308 << endl;
2309 }
2310/*
2311 const label uDegree(uBasis_.degree());
2312 const label vDegree(vBasis_.degree());
2313 const label uNCPs(uBasis_.nCPs());
2314 const label vNCPs(vBasis_.nCPs());
2315
2316 OFstream surfaceFileUBases(dirName/fileName+"UBases");
2317 OFstream surfaceFileVBases(dirName/fileName+"VBases");
2318
2319 forAll(*this, ptI)
2320 {
2321 const scalar& u(u_[ptI]);
2322 const scalar& v(v_[ptI]);
2323 scalarField uBasesValues(uNCPs);
2324 scalarField vBasesValues(vNCPs);
2325
2326 surfaceFileUBases << u << " ";
2327 surfaceFileVBases << v << " ";
2328
2329 forAll(CPs_, CPI)
2330 {
2331 basesValues[CPI] = basis_.basisValue(CPI, degree, u);
2332 curveFileBases << basesValues[CPI] << " ";
2333 }
2334
2335 curveFileBases << endl;
2336 }
2337*/
2338 }
2339}
2340
2341
2343(
2344 const fileName vtkDirName,
2345 const fileName vtkFileName
2346)
2347{
2348 if (Pstream::master())
2349 {
2350 if (vtkFileName.has_ext())
2351 {
2353 << "Do not supply a file extension."
2354 << exit(FatalError);
2355 }
2356
2357 // Build the surface.
2358 buildSurface();
2359
2360 // Write the surface as vtk.
2361 OFstream surfaceFile(vtkFileName);
2362 pointField& surfacePoints(*this);
2363 faceList surfaceFaces((nUPts_ - 1)*(nUPts_ - 1), face(4));
2364
2365 for (label fuI = 0; fuI < (nUPts_ - 1); fuI++)
2366 {
2367 for (label fvI = 0; fvI < (nVPts_ - 1); fvI++)
2368 {
2369 const label fI(fuI*(nUPts_ - 1) + fvI);
2370 face& surfaceFace(surfaceFaces[fI]);
2371
2372 surfaceFace[0] = (fuI)*nVPts_ + (fvI);
2373 surfaceFace[1] = (fuI + 1)*nVPts_ + (fvI);
2374 surfaceFace[2] = (fuI + 1)*nVPts_ + (fvI + 1);
2375 surfaceFace[3] = (fuI)*nVPts_ + (fvI + 1);
2376 }
2377 }
2378
2379 surfaceWriters::vtkWriter writer;
2380
2381 writer.open
2382 (
2383 surfacePoints,
2384 surfaceFaces,
2385 vtkDirName/vtkFileName,
2386 false
2387 );
2388
2389 writer.close();
2390
2391 // Write the control surface as vtk.
2392 fileName vtkCPFileName(vtkFileName+"CPs");
2393 pointField surfaceCPPoints(CPs_);
2394 const label uNCPs(uBasis_.nCPs());
2395 const label vNCPs(vBasis_.nCPs());
2396 faceList surfaceCPFaces((uNCPs-1)*(vNCPs-1), face(4));
2397
2398 for (label fvCPI=0; fvCPI<(vNCPs-1); fvCPI++)
2399 {
2400 for (label fuCPI=0; fuCPI<(uNCPs-1); fuCPI++)
2401 {
2402 const label fCPI(fvCPI*(uNCPs-1) + fuCPI);
2403 face& surfaceCPFace(surfaceCPFaces[fCPI]);
2404
2405 surfaceCPFace[0] = (fvCPI)*uNCPs + (fuCPI);
2406 surfaceCPFace[1] = (fvCPI + 1)*uNCPs + (fuCPI);
2407 surfaceCPFace[2] = (fvCPI + 1)*uNCPs + (fuCPI + 1);
2408 surfaceCPFace[3] = (fvCPI)*uNCPs + (fuCPI + 1);
2409 }
2410 }
2411
2412 writer.open
2413 (
2414 surfaceCPPoints,
2415 surfaceCPFaces,
2416 vtkDirName/vtkCPFileName,
2417 false
2418 );
2419
2420 writer.close();
2421 }
2422}
2423
2424
2425// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2426
2427} // End namespace Foam
2428
2429// ************************************************************************* //
scalar delta
#define R(A, B, C, D, E, F, K, M)
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
scalar surfaceDerivativeCP(const scalar u, const scalar v, const label CPI) const
Surface derivative wrt the weight of CPI at point u,v.
List< scalarList > genEquidistant(const label nUPts=100, const label nVPts=100, const label lenAcc=25, const label maxIter=10, const label spacingCorrInterval=-1, const scalar tolerance=1.e-5)
Generate points on the surface which are evenly spaced in cartesian.
const labelList & getBoundaryCPIs()
vector surfaceDerivativeUV(const scalar u, const scalar v) const
Surface second derivative wrt u and v at point u,v.
bool checkRangeU(const scalar u, const label CPI, const label uDegree) const
void writeVTK(const fileName vtkDirName, const fileName vtkFileName)
void setCPs(const List< vector > &CPs)
scalar lengthU(const label vIConst, const label uIStart, const label uIEnd) const
void makeEquidistant(const label lenAcc=25, const label maxIter=10, const label spacingCorrInterval=-1, const scalar tolerance=1.e-5)
scalarList findClosestSurfacePoint(const vector &targetPoint, const label maxIter=100, const scalar tolerance=1.e-6)
vector surfacePoint(const scalar &u, const scalar &v)
vector surfaceDerivativeV(const scalar u, const scalar v) const
Surface derivative wrt v at point u,v.
void setName(const word &name)
scalar lengthDerivativeV(const scalar uConst, const scalar vStart, const scalar vEnd, const label nPts) const
Surface derivative wrt v length along u=const contour range.
void setNrmOrientation(const vector &givenNrm, const scalar u, const scalar v)
scalar lengthV(const label uIConst, const label vIStart, const label vIEnd) const
NURBS3DSurface(const List< vector > &CPs, const label nPointsU, const label nPointsV, const NURBSbasis &uBasis, const NURBSbasis &vBasis, const word name="NURBS3DSurface")
Construct from number of control points and basis functions.
const labelList & getBoundaryCPIDs()
Get IDs of boundary control points.
void write()
Write curve to file.
vector surfaceDerivativeUU(const scalar u, const scalar v) const
Surface second derivative wrt u at point u,v.
vector surfaceDerivativeVV(const scalar u, const scalar v) const
Surface second derivative wrt v at point u,v.
const vector nrm(scalar u, scalar v)
bool checkRangeV(const scalar v, const label CPI, const label vDegree) const
vector surfaceDerivativeW(const scalar u, const scalar v, const label CPI) const
Surface derivative wrt WI at point u,v.
vector surfaceDerivativeU(const scalar u, const scalar v) const
Surface derivative wrt u at point u,v.
scalar lengthDerivativeU(const scalar vConst, const scalar uStart, const scalar uEnd, const label nPts) const
Surface derivative wrt u length along v=const contour range.
void setWeights(const scalarList &weights)
void flipNrmOrientation()
Flip the orientation of the nrm.
bool checkRangeUV(const scalar v, const scalar u, const label CPI, const label uDegree, const label vDegree) const
const label & whichBoundaryCPI(const label &globalCPI)
Get the boundary CP ID based on the global CP ID.
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
Definition NURBSbasis.H:52
scalar basisDerivativeUU(const label iCP, const label degree, const scalar u) const
Basis second derivative w.r.t u.
Definition NURBSbasis.C:224
scalar basisDerivativeU(const label iCP, const label degree, const scalar u) const
Basis derivative w.r.t u.
Definition NURBSbasis.C:183
const label & nCPs() const
Definition NURBSbasisI.H:39
const label & degree() const
Definition NURBSbasisI.H:33
scalar basisValue(const label iCP, const label degree, const scalar u) const
Basis value.
Definition NURBSbasis.C:133
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
void size(const label n)
Definition UList.H:118
vector & operator[](const label i)
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static const Form zero
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
bool has_ext() const
Various checks for extensions.
Definition stringI.H:43
A surfaceWriter for VTK legacy (.vtk) or XML (.vtp) format.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
rDeltaTY field()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define WarningInFunction
Report a warning using Foam::Warning.
const expr V(m.psi().mesh().V())
const dimensionedScalar c
Speed of light in a vacuum.
const wordList surface
Standard surface field types (scalar, vector, tensor, etc).
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299