Go to the documentation of this file.
39 #ifndef NURBS3DSurface_H
40 #define NURBS3DSurface_H
100 label nrmOrientation_;
109 label sgn(
const scalar val)
const;
111 scalar abs(
const scalar val)
const;
113 label mod(
const label
x,
const label interval)
const;
115 void setCPUVLinking();
131 const scalar minVal=1
e-7,
132 const scalar maxVal=0.999999
138 const scalar minVal=1
e-7,
139 const scalar maxVal=0.999999
152 const label spacingCorrInterval,
153 const scalar tolerance
164 const label nPointsU,
165 const label nPointsV,
176 const label nPointsU,
177 const label nPointsV,
187 const label nPointsU,
188 const label nPointsV,
201 const label nPointsU,
202 const label nPointsV,
218 const label nPointsU,
219 const label nPointsV,
233 const label nPointsU,
234 const label nPointsV,
282 const label lenAcc=25,
283 const label maxIter=10,
284 const label spacingCorrInterval=-1,
285 const scalar tolerance=1.
e-5
294 const vector& targetPoint,
295 const label maxIter=100,
296 const scalar tolerance=1.
e-6
302 const label maxIter=100,
303 const scalar tolerance=1.
e-6
308 const vector& targetPoint,
309 const scalar& uInitGuess,
310 const scalar& vInitGuess,
311 const label maxIter=100,
312 const scalar tolerance=1.
e-6
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
338 bool checkRangeU(
const scalar u,
const label CPI)
const;
347 bool checkRangeV(
const scalar v,
const label CPI)
const;
380 scalar
lengthU(
const label vIConst)
const;
382 scalar
lengthU(
const scalar vConst)
const;
399 scalar
lengthV(
const label uIConst)
const;
401 scalar
lengthV(
const scalar uConst)
const;
518 inline label nrmOrientation()
const
520 return nrmOrientation_;
526 return givenInitNrm_;
546 void write(
const word fileName);
548 void write(
const fileName dirName,
const fileName fileName);
554 void writeWParses(
const fileName dirName,
const fileName fileName);
556 void writeVTK(
const fileName vtkDirName,
const fileName vtkFileName);
List< label > labelList
A List of labels.
vector surfacePoint(const scalar &u, const scalar &v)
~NURBS3DSurface()=default
const labelList & getCPsUCPIs() const
const label & whichBoundaryCPI(const label &globalCPI)
vector surfaceDerivativeU(const scalar u, const scalar v) const
vector surfaceDerivativeUU(const scalar u, const scalar v) const
A class for handling words, derived from Foam::string.
scalar lengthDerivativeV(const scalar uConst, const scalar vStart, const scalar vEnd, const label nPts) const
A class for handling file names.
A class for managing temporary objects.
void flipNrmOrientation()
const word & getName() const
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
NURBS3DSurface(const List< vector > &CPs, const label nPointsU, const label nPointsV, const NURBSbasis &uBasis, const NURBSbasis &vBasis, const word name="NURBS3DSurface")
void setName(const word &name)
void writeVTK(const fileName vtkDirName, const fileName vtkFileName)
const scalarList & getWeights() const
scalar lengthDerivativeU(const scalar vConst, const scalar uStart, const scalar uEnd, const label nPts) const
bool checkRangeV(const scalar v, const label CPI, const label vDegree) const
vector surfaceDerivativeW(const scalar u, const scalar v, const label CPI) const
bool checkRangeUV(const scalar v, const scalar u, const label CPI, const label uDegree, const label vDegree) const
const vector & givenInitNrm() const
bool checkRangeU(const scalar u, const label CPI, const label uDegree) const
#define R(A, B, C, D, E, F, K, M)
Generic templated field type.
const scalarList & getParametricCoordinatesU() const
vector surfaceDerivativeV(const scalar u, const scalar v) const
scalarList findClosestSurfacePoint(const vector &targetPoint, const label maxIter=100, const scalar tolerance=1.e-6)
void setNrmOrientation(const vector &givenNrm, const scalar u, const scalar v)
scalar lengthV(const label uIConst, const label vIStart, const label vIEnd) const
void setWeights(const scalarList &weights)
const labelList & getBoundaryCPIs()
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
vector surfaceDerivativeUV(const scalar u, const scalar v) const
const vector nrm(scalar u, scalar v)
scalar lengthU(const label vIConst, const label uIStart, const label uIEnd) const
const labelList & getCPsVCPIs() const
const dimensionedScalar e
vector surfaceDerivativeVV(const scalar u, const scalar v) const
word name(const expressions::valueTypeCode typeCode)
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
const scalarList & getParametricCoordinatesV() const
const labelList & getBoundaryCPIDs()
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)
const NURBSbasis & getBasisFunctionU() const
scalar surfaceDerivativeCP(const scalar u, const scalar v, const label CPI) const
void setCPs(const List< vector > &CPs)
label nrmOrientation() const