NURBS3DSurface.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 -------------------------------------------------------------------------------
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 Class
29  Foam::NURBS3DSurface
30 
31 Description
32  A NURBS 3D surface
33 
34 SourceFiles
35  NURBS3DSurface.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef NURBS3DSurface_H
40 #define NURBS3DSurface_H
41 
42 #include "NURBSbasis.H"
43 #include "vectorField.H"
44 #include "vector2DField.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class NURBS3DSurface Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class NURBS3DSurface
56 :
57  public vectorField
58 {
59  // Private Data
60 
61  enum paramType
62  {
63  PARAMU=0,
64  PARAMV,
65  NPARAMS
66  };
67 
68  enum nrmOrientation
69  {
70  ALIGNED=1,
71  OPPOSED=-1
72  };
73 
74  //- Held in sets of const v.
75  List<vector> CPs_;
76 
77  //- Held in sets of const u.
78  scalarList u_;
79 
80  scalarList v_;
81 
82  scalarList weights_;
83 
84  label nUPts_;
85 
86  label nVPts_;
87 
88  word name_;
89 
90  NURBSbasis uBasis_;
91 
92  NURBSbasis vBasis_;
93 
94  vector givenInitNrm_;
95 
96  labelList CPsUCPIs_;
97 
98  labelList CPsVCPIs_;
99 
100  label nrmOrientation_;
101 
102  autoPtr<labelList> boundaryCPIDs_;
103 
104  autoPtr<labelList> whichBoundaryCPID_;
105 
106 
107  // Private Member Functions
108 
109  label sgn(const scalar val) const;
110 
111  scalar abs(const scalar val) const;
112 
113  label mod(const label x, const label interval) const;
114 
115  void setCPUVLinking();
116 
117  void setUniformUV
118  (
119  scalarList& u,
120  scalarList& v,
121  const label nUPts,
122  const label nVPts
123  ) const;
124 
125  void setUniformUV();
126 
127  bool bound
128  (
129  scalar& u,
130  scalar& v,
131  const scalar minVal=1e-7,
132  const scalar maxVal=0.999999
133  ) const;
134 
135  bool boundDirection
136  (
137  scalar& u,
138  const scalar minVal=1e-7,
139  const scalar maxVal=0.999999
140  ) const;
141 
142  // Function structure assumes that range specified for param R
143  // should range from 0->1 and is sized correctly for nPts.
144  // Pts are set in equidistant manner along fixed value of other param S.
145  void setEquidistantR
146  (
147  scalarList& R,
148  const scalar SHeld,
149  const label paramR,
150  const label lenAcc,
151  const label maxIter,
152  const label spacingCorrInterval,
153  const scalar tolerance
154  ) const;
155 
156 public:
157 
158  // Constructors
159 
160  //- Construct from number of control points and basis functions
162  (
163  const List<vector>& CPs,
164  const label nPointsU,
165  const label nPointsV,
166  const NURBSbasis& uBasis,
167  const NURBSbasis& vBasis,
168  const word name="NURBS3DSurface"
169  );
170 
171  //- Construct from number of control points, weights and basis functions
173  (
174  const List<vector>& CPs,
175  const List<scalar>& weights,
176  const label nPointsU,
177  const label nPointsV,
178  const NURBSbasis& uBasis,
179  const NURBSbasis& vBasis,
180  const word name="NURBS3DSurface"
181  );
182 
183  //- Construct from control points, basis degree and number of points
185  (
186  const List<vector>& CPs,
187  const label nPointsU,
188  const label nPointsV,
189  const label uDegree,
190  const label vDegree,
191  const label nCPsU,
192  const label nCPsV,
193  const word name="NURBS3DSurface"
194  );
195 
196  //- Construct from control points, basis degree, knots and number of
197  // points
199  (
200  const List<vector>& CPs,
201  const label nPointsU,
202  const label nPointsV,
203  const label uDegree,
204  const label vDegree,
205  const label nCPsU,
206  const label nCPsV,
207  const scalarField& knotsU,
208  const scalarField& knotsV,
209  const word name="NURBS3DSurface"
210  );
211 
212  //- Construct from control points, weights, basis degree, and number of
213  // points
215  (
216  const List<vector>& CPs,
217  const List<scalar>& weights,
218  const label nPointsU,
219  const label nPointsV,
220  const label uDegree,
221  const label vDegree,
222  const label nCPsU,
223  const label nCPsV,
224  const word name="NURBS3DSurface"
225  );
226 
227  //- Construct from control points, weights, basis degree, knots and
228  // number of points
230  (
231  const List<vector>& CPs,
232  const List<scalar>& weights,
233  const label nPointsU,
234  const label nPointsV,
235  const label uDegree,
236  const label vDegree,
237  const label nCPsU,
238  const label nCPsV,
239  const scalarField& knotsU,
240  const scalarField& knotsV,
241  const word name="NURBS3DSurface"
242  );
243 
244  //- Construct as copy
246 
247 
248  //- Destructor
249 
250  ~NURBS3DSurface() = default;
251 
252  // Member Functions
253 
254  // Set Functions
255 
256  void setNrmOrientation
257  (
258  const vector& givenNrm,
259  const scalar u,
260  const scalar v
261  );
262 
263  //- Flip the orientation of the nrm.
264  void flipNrmOrientation();
265 
266  void setCPs(const List<vector>& CPs);
267 
268  void setWeights(const scalarList& weights);
269 
270  void setName(const word& name);
271 
272  void buildSurface();
273 
274  void invertU();
275 
276  void invertV();
277 
278  void invertUV();
279 
280  void makeEquidistant
281  (
282  const label lenAcc=25,
283  const label maxIter=10,
284  const label spacingCorrInterval=-1,
285  const scalar tolerance=1.e-5
286  );
287 
288  // Calculation Functions
289 
290  vector surfacePoint( const scalar& u, const scalar& v);
291 
293  (
294  const vector& targetPoint,
295  const label maxIter=100,
296  const scalar tolerance=1.e-6
297  );
298 
300  (
301  const vectorField& targetPoints,
302  const label maxIter=100,
303  const scalar tolerance=1.e-6
304  );
305 
307  (
308  const vector& targetPoint,
309  const scalar& uInitGuess,
310  const scalar& vInitGuess,
311  const label maxIter=100,
312  const scalar tolerance=1.e-6
313  );
314 
315  const vector nrm(scalar u, scalar v);
316 
317  //- Generate points on the surface which are evenly spaced in cartesian
318  // coordinate distances.
320  (
321  const label nUPts=100,
322  const label nVPts=100,
323  const label lenAcc=25,
324  const label maxIter=10,
325  const label spacingCorrInterval=-1,
326  const scalar tolerance=1.e-5
327  );
328 
329  // Location Functions
330 
331  bool checkRangeU
332  (
333  const scalar u,
334  const label CPI,
335  const label uDegree
336  ) const;
337 
338  bool checkRangeU(const scalar u, const label CPI) const;
339 
340  bool checkRangeV
341  (
342  const scalar v,
343  const label CPI,
344  const label vDegree
345  ) const;
346 
347  bool checkRangeV(const scalar v, const label CPI) const;
348 
349  bool checkRangeUV
350  (
351  const scalar v,
352  const scalar u,
353  const label CPI,
354  const label uDegree,
355  const label vDegree
356  ) const;
357 
358  bool checkRangeUV
359  (
360  const scalar v,
361  const scalar u,
362  const label CPI
363  ) const;
364 
365  scalar lengthU
366  (
367  const label vIConst,
368  const label uIStart,
369  const label uIEnd
370  ) const;
371 
372  scalar lengthU
373  (
374  const scalar vConst,
375  const scalar uStart,
376  const scalar uEnd,
377  const label nPts
378  ) const;
379 
380  scalar lengthU(const label vIConst) const;
381 
382  scalar lengthU(const scalar vConst) const;
383 
384  scalar lengthV
385  (
386  const label uIConst,
387  const label vIStart,
388  const label vIEnd
389  ) const;
390 
391  scalar lengthV
392  (
393  const scalar uConst,
394  const scalar vStart,
395  const scalar vEnd,
396  const label nPts
397  ) const;
398 
399  scalar lengthV(const label uIConst) const;
400 
401  scalar lengthV(const scalar uConst) const;
402 
403  // Derivative Functions
404 
405  //- Surface derivative wrt u at point u,v.
406  vector surfaceDerivativeU(const scalar u, const scalar v) const;
407 
408  //- Surface derivative wrt v at point u,v.
409  vector surfaceDerivativeV(const scalar u, const scalar v) const;
410 
411  //- Surface second derivative wrt u and v at point u,v.
412  vector surfaceDerivativeUV(const scalar u, const scalar v) const;
413 
414  //- Surface second derivative wrt u at point u,v.
415  vector surfaceDerivativeUU(const scalar u, const scalar v) const;
416 
417  //- Surface second derivative wrt v at point u,v.
418  vector surfaceDerivativeVV(const scalar u, const scalar v) const;
419 
420  //- Surface derivative wrt the weight of CPI at point u,v.
421  scalar surfaceDerivativeCP
422  (
423  const scalar u,
424  const scalar v,
425  const label CPI
426  ) const;
427 
428  //- Surface derivative wrt WI at point u,v.
430  (
431  const scalar u,
432  const scalar v,
433  const label CPI
434  ) const;
435 
436  //- Surface derivative wrt u length along v=const contour range.
437  scalar lengthDerivativeU
438  (
439  const scalar vConst,
440  const scalar uStart,
441  const scalar uEnd,
442  const label nPts
443  ) const;
444 
445  //- Surface derivative wrt v length along u=const contour range.
446  scalar lengthDerivativeV
447  (
448  const scalar uConst,
449  const scalar vStart,
450  const scalar vEnd,
451  const label nPts
452  ) const;
453 
454 
455  // Access Functions
456 
457  //- Get basis function.
458  inline const NURBSbasis& getBasisFunctionU() const
459  {
460  return uBasis_;
461  }
462 
463  inline const NURBSbasis& getBasisFunctionV() const
464  {
465  return vBasis_;
466  }
467 
468  //- Get CPs.
469  inline const List<vector>& getCPs() const
470  {
471  return CPs_;
472  }
473 
474  //- Get weights.
475  inline const scalarList& getWeights() const
476  {
477  return weights_;
478  }
479 
480  //- Get parametric coordinates.
481  inline const scalarList& getParametricCoordinatesU() const
482  {
483  return u_;
484  }
485 
486  inline const scalarList& getParametricCoordinatesV() const
487  {
488  return v_;
489  }
490 
491  //- Get name.
492  inline const word& getName() const
493  {
494  return name_;
495  }
496 
497  //- Get number of point in u direction
498  inline label getNPtsU() const
499  {
500  return nUPts_;
501  }
502 
503  //- Get number of point in u direction
504  inline label getNPtsV() const
505  {
506  return nVPts_;
507  }
508 
509  //- Get IDs of boundary control points
510  const labelList& getBoundaryCPIDs();
511 
512  const labelList& getBoundaryCPIs();
513 
514  //- Get the boundary CP ID based on the global CP ID
515  const label& whichBoundaryCPI(const label& globalCPI);
516 
517  //- Return the nrm sgn relation to the u=0 nrm.
518  inline label nrmOrientation() const
519  {
520  return nrmOrientation_;
521  }
523  //- Return the initial nrmal given to compare to the Curve's nrmals.
524  inline const vector& givenInitNrm() const
525  {
526  return givenInitNrm_;
527  }
528 
529  //- Return ID in u direction for a given cp ID
530  const labelList& getCPsUCPIs() const
531  {
532  return CPsUCPIs_;
533  }
534 
535  //- Return ID in v direction for a given cp ID
536  const labelList& getCPsVCPIs() const
537  {
538  return CPsVCPIs_;
539  }
540 
541  // Write Functions
542 
543  //- Write curve to file
544  void write();
545 
546  void write(const word fileName);
547 
548  void write(const fileName dirName, const fileName fileName);
549 
550  void writeWParses();
552  void writeWParses(const word fileName);
553 
554  void writeWParses(const fileName dirName, const fileName fileName);
555 
556  void writeVTK(const fileName vtkDirName, const fileName vtkFileName);
557 };
558 
559 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
560 
561 } // End namespace Foam
562 
563 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
564 
565 #endif
566 
567 // ************************************************************************* //
scalar lengthU(const label vIConst, const label uIStart, const label uIEnd) const
const labelList & getCPsVCPIs() const
Return ID in v direction for a given cp ID.
label getNPtsU() const
Get number of point in u direction.
vector surfaceDerivativeUV(const scalar u, const scalar v) const
Surface second derivative wrt u and v at point u,v.
A class for handling file names.
Definition: fileName.H:72
vector surfaceDerivativeVV(const scalar u, const scalar v) const
Surface second derivative wrt v at point u,v.
const scalarList & getParametricCoordinatesU() const
Get parametric coordinates.
const NURBSbasis & getBasisFunctionV() const
void makeEquidistant(const label lenAcc=25, const label maxIter=10, const label spacingCorrInterval=-1, const scalar tolerance=1.e-5)
const List< vector > & getCPs() const
Get CPs.
const labelList & getBoundaryCPIDs()
Get IDs of boundary control points.
scalar lengthV(const label uIConst, const label vIStart, const label vIEnd) const
scalar surfaceDerivativeCP(const scalar u, const scalar v, const label CPI) const
Surface derivative wrt the weight of CPI at point u,v.
A NURBS 3D surface.
void setCPs(const List< vector > &CPs)
vector surfacePoint(const scalar &u, const scalar &v)
void write()
Write curve to file.
~NURBS3DSurface()=default
Destructor.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
const vector nrm(scalar u, scalar v)
const word & getName() const
Get name.
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
Definition: NURBSbasis.H:51
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 scalarList & getWeights() const
Get weights.
const scalarList & getParametricCoordinatesV() const
vector surfaceDerivativeW(const scalar u, const scalar v, const label CPI) const
Surface derivative wrt WI at point u,v.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const vector & givenInitNrm() const
Return the initial nrmal given to compare to the Curve&#39;s nrmals.
bool checkRangeUV(const scalar v, const scalar u, const label CPI, const label uDegree, const label vDegree) const
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.
label getNPtsV() const
Get number of point in u direction.
bool checkRangeU(const scalar u, const label CPI, const label uDegree) const
const NURBSbasis & getBasisFunctionU() const
Get basis function.
const labelList & getCPsUCPIs() const
Return ID in u direction for a given cp ID.
bool checkRangeV(const scalar v, const label CPI, const label vDegree) const
label nrmOrientation() const
Return the nrm sgn relation to the u=0 nrm.
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.
scalarList findClosestSurfacePoint(const vector &targetPoint, const label maxIter=100, const scalar tolerance=1.e-6)
#define R(A, B, C, D, E, F, K, M)
vector surfaceDerivativeU(const scalar u, const scalar v) const
Surface derivative wrt u at point u,v.
vector surfaceDerivativeUU(const scalar u, const scalar v) const
Surface second derivative wrt u at point u,v.
const label & whichBoundaryCPI(const label &globalCPI)
Get the boundary CP ID based on the global CP ID.
void writeVTK(const fileName vtkDirName, const fileName vtkFileName)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
void setWeights(const scalarList &weights)
vector surfaceDerivativeV(const scalar u, const scalar v) const
Surface derivative wrt v at point u,v.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Foam::vector2DField.
void setNrmOrientation(const vector &givenNrm, const scalar u, const scalar v)
const labelList & getBoundaryCPIs()
void flipNrmOrientation()
Flip the orientation of the nrm.
void setName(const word &name)
Namespace for OpenFOAM.
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.