NURBS3DCurve.H
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-------------------------------------------------------------------------------
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
29Class
30 Foam::NURBS3DCurve
31
32Description
33 A NURBS 3D curve
34
35SourceFiles
36 NURBS3DCurve.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef NURBS3DCurve_H
41#define NURBS3DCurve_H
42
43#include "NURBSbasis.H"
44#include "vectorField.H"
45#include "OFstream.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52/*---------------------------------------------------------------------------*\
53 Class NURBS3DCurve Declaration
54\*---------------------------------------------------------------------------*/
56class NURBS3DCurve
57:
58 public vectorField
59{
60 // Private Data
61
62 enum nrmOrientation
63 {
64 ALIGNED = 1,
65 OPPOSED = -1
66 };
67
68 //- NURBS basis function
69 List<vector> CPs_;
70 List<scalar> weights_;
71 scalarList u_;
72 word name_;
73 const NURBSbasis& basis_;
74
75 vector givenInitNrm_;
76
77 label nrmOrientation_;
78
79private:
80
81 // Private Member Functions
82
83 label sgn(const scalar val) const;
84
85 scalar abs(const scalar val) const;
86
87 label mod(const label x, const label interval) const;
88
89 void setUniformU();
90
91 // Bound parametric coor to certain limits.
92 bool bound
93 (
94 scalar& u,
95 const scalar minVal,
96 const scalar maxVal
97 ) const;
98
99 // Function structure assumes that U is from 0->1 and pre-sized to nPts.
100 void setEquidistantU
101 (
102 scalarList& U,
103 const label lenAcc,
104 const label maxIter,
105 const label spacingCorrInterval,
106 const scalar tolerance
107 ) const;
108
109
110public:
111
112 // Constructors
113
114 //- Construct from control points, weights and basis order and
115 //- parametric coordinates
117 (
118 const NURBSbasis& basis,
119 const List<vector>& CPs,
120 const List<scalar>& weights,
121 const scalarField& u,
122 const label nPts,
123 const word name="NURBS3DCurve"
124 );
125
126 //- Construct from control points, weights and basis order.
127 // Uniform coordinate distribution is implied
129 (
130 const NURBSbasis& basis,
131 const List<vector>& CPs,
132 const List<scalar>& weights,
133 const label nPts,
134 const word name="NURBS3DCurve"
135 );
136
137 //- Construct from control points and basis order
138 // Uniform coordinate distribution and unit weights are implied
140 (
141 const NURBSbasis& basis,
142 const List<vector>& CPs,
143 const label nPts,
144 const word name="NURBS3DCurve"
145 );
146
147 //- Construct as copy
149
150
151 //- Destructor
152 ~NURBS3DCurve() = default;
153
154
155 // Member Functions
156
157 // Set Functions
158
159 //- Take a given normal and use to determine if
160 //- NURBS normals should be reversed.
161 //- Computation taken from u = 0.
163 (
164 const vector& givenNrm,
165 const vector& givenTan
166 );
167
169 (
170 const vector& givenNrm,
171 const scalar zVal
172 );
173
174 //- Flip the orientation of the nrm.
175 void flipNrmOrientation();
176
177 //- Set CPs.
178 void setCPs(const List<vector>& CPs);
179
180 //- Set weights.
181 void setWeights(const List<scalar>& weights);
182
183 //- Set name.
184 void setName(const word& name);
185
186 //- Build curve.
187 void buildCurve();
188
189 //- Invert CPs order and re-build curve.
190 void invert();
191
192 //- Insert a knot by re-computing the control points.
193 // The basis' insertKnot function must have beem called first to
194 // correctly provide the required info.
195 void insertKnot
196 (
197 const scalarField& oldKnots,
198 const scalar uBar,
199 const label kInsert
200 );
201
202 //- Make curve points equidistant in cartesian space.
203 void makeEquidistant
204 (
205 const label lenAcc = 25,
206 const label maxIter = 10,
207 const label spacingCorrInterval = -1,
208 const scalar tolerance = 1.e-5
209 );
210
211
212 // Point Calc Functions
213
214 //- Curve point cartesian coordinates at ui.
215 vector curvePoint(const scalar u) const;
216
217 //- Find curve point which is closest to given point
218 //- using Newton-Raphson. Returns the param coordinate.
220 (
221 const vector& targetPoint,
222 const label maxIter = 1000,
223 const scalar tolerance = 1.e-13
224 );
225
226 //- Find curve point which is closest to given point
227 //- using Newton-Raphson. Returns the param coordinate.
229 (
230 const vector& targetPoint,
231 const scalar initGuess,
232 const label maxIter = 1000,
233 const scalar tolerance = 1.e-13
234 );
235
236 //- Find the normal to the curve, with the option of forcing a z-plane.
237 const vector nrm3D(const vector& refTan, const scalar u) const;
238
239 const vector nrm2D(const scalar zVal, const scalar u) const;
240
241 //- Generate points on the curve which are evenly spaced in cartesian
242 //- coordinate distances.
244 (
245 const label nPts = 100,
246 const label lenAcc = 25,
247 const label maxIter = 10,
248 const label spacingCorrInterval=-1,
249 const scalar tolerance = 1.e-5
250 );
251
252
253 // Location Functions
254
255 //- Check if given parametric coordinate u and CP are linked.
256 bool checkRange
257 (
258 const scalar u,
259 const label CPI,
260 const label degree
261 ) const;
262
263 bool checkRange
264 (
265 const scalar u,
266 const label CPI
267 ) const;
268
269 //- Calculate Length from starting to ending indices via
270 //- computational evaluation using trapezoid rule.
271 scalar length
272 (
273 const label uIStart,
274 const label uIEnd
275 ) const;
276
277 //- Calculate length from starting to ending parametric coordinates via
278 //- computational evaluation using trapezoid rule.
279 scalar length
280 (
281 const scalar uStart,
282 const scalar uEnd,
283 const label nPts
284 ) const;
285
286 //- Calculate length for the entire curve length.
287 scalar length() const;
288
289
290 // Derivative Functions
291
292 //- Curve derivative wrt u at point ui.
293 vector curveDerivativeU(const scalar u) const;
294
295 //- Curve second derivative wrt u at point ui.
296 vector curveDerivativeUU(const scalar u) const;
297
298 //- Curve derivative wrt b at point ui:
299 //- scalar since dx/dX = dy/dY = dz/dZ.
300 scalar curveDerivativeCP
301 (
302 const scalar u,
303 const label CPI
304 );
305
306 //- Curve derivative wrt CPII at point u.
308 (
309 const scalar u,
310 const label CPI
311 );
312
313 //- Calculate length from starting to ending indices via
314 //- computational evaluation using trapezoid rule
315 scalar lengthDerivativeU
316 (
317 const scalar uStart,
318 const scalar uEnd,
319 const label nPts
320 );
321
322
323 // Access Functions
324
325 //- Get basis function.
326 inline const NURBSbasis& getBasisFunction() const
327 {
328 return basis_;
329 }
330
331 //- Get CPs.
332 inline const List<vector>& getCPs() const
333 {
334 return CPs_;
335 }
336
337 //- Get weights.
338 inline const List<scalar>& getWeights() const
339 {
340 return weights_;
341 }
342
343 //- Get parametric coordinates.
344 inline const List<scalar>& getParametricCoordinates() const
345 {
346 return u_;
347 }
348
349 //- Get name.
350 inline const word& getName() const
351 {
352 return name_;
353 }
354
355 //- Return the nrm sgn relation to the u=0 nrm.
356 inline label nrmOrientation() const
357 {
358 return nrmOrientation_;
359 }
360
361 //- Return the initial normal given to compare to the Curve's normals.
362 inline const vector& givenInitNrm() const
363 {
364 return givenInitNrm_;
365 }
366
367
368 // Write Functions
369
370 //- Write curve to file.
371 void write();
372
373 void write
374 (
375 const word fileName
376 );
377
378 void write
379 (
380 const fileName dirName,
381 const fileName fileName
382 );
383
384 void writeWParses();
385
386 void writeWParses
387 (
388 const word fileName
389 );
390
391 void writeWParses
392 (
393 const fileName dirName,
394 const fileName fileName
395 );
396};
397
398
399// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
400
401} // End namespace Foam
402
403// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404
405#endif
406
407// ************************************************************************* //
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A NURBS 3D curve.
Definition: NURBS3DCurve.H:58
~NURBS3DCurve()=default
Destructor.
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
void insertKnot(const scalarField &oldKnots, const scalar uBar, const label kInsert)
Insert a knot by re-computing the control points.
Definition: NURBS3DCurve.C:646
void setNrm3DOrientation(const vector &givenNrm, const vector &givenTan)
Definition: NURBS3DCurve.C:291
vector curveDerivativeUU(const scalar u) const
Curve second derivative wrt u at point ui.
Definition: NURBS3DCurve.C:825
NURBS3DCurve(const NURBS3DCurve &)
Construct as copy.
void setCPs(const List< vector > &CPs)
Set CPs.
Definition: NURBS3DCurve.C:359
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
vector curveDerivativeWeight(const scalar u, const label CPI)
Curve derivative wrt CPII at point u.
Definition: NURBS3DCurve.C:886
scalar length() const
Calculate length for the entire curve length.
Definition: NURBS3DCurve.C:792
void setName(const word &name)
Set name.
Definition: NURBS3DCurve.C:371
const List< scalar > & getWeights() const
Get weights.
Definition: NURBS3DCurve.H:337
void invert()
Invert CPs order and re-build curve.
Definition: NURBS3DCurve.C:408
const List< scalar > & getParametricCoordinates() const
Get parametric coordinates.
Definition: NURBS3DCurve.H:343
scalar findClosestCurvePoint(const vector &targetPoint, const label maxIter=1000, const scalar tolerance=1.e-13)
Definition: NURBS3DCurve.C:484
label nrmOrientation() const
Return the nrm sgn relation to the u=0 nrm.
Definition: NURBS3DCurve.H:355
void buildCurve()
Build curve.
Definition: NURBS3DCurve.C:377
const NURBSbasis & getBasisFunction() const
Get basis function.
Definition: NURBS3DCurve.H:325
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
const vector nrm2D(const scalar zVal, const scalar u) const
Definition: NURBS3DCurve.C:631
scalar lengthDerivativeU(const scalar uStart, const scalar uEnd, const label nPts)
Definition: NURBS3DCurve.C:912
const vector & givenInitNrm() const
Return the initial normal given to compare to the Curve's normals.
Definition: NURBS3DCurve.H:361
void write()
Write curve to file.
Definition: NURBS3DCurve.C:954
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
const word & getName() const
Get name.
Definition: NURBS3DCurve.H:349
vector curveDerivativeU(const scalar u) const
Curve derivative wrt u at point ui.
Definition: NURBS3DCurve.C:800
void setNrm2DOrientation(const vector &givenNrm, const scalar zVal)
Definition: NURBS3DCurve.C:317
void flipNrmOrientation()
Flip the orientation of the nrm.
Definition: NURBS3DCurve.C:346
void setWeights(const List< scalar > &weights)
Set weights.
Definition: NURBS3DCurve.C:365
vector curvePoint(const scalar u) const
Curve point cartesian coordinates at ui.
Definition: NURBS3DCurve.C:457
scalar curveDerivativeCP(const scalar u, const label CPI)
Definition: NURBS3DCurve.C:864
const List< vector > & getCPs() const
Get CPs.
Definition: NURBS3DCurve.H:331
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
Definition: NURBSbasis.H:56
A class for handling file names.
Definition: fileName.H:76
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59