NURBS3DVolume Class Referenceabstract

NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing support for different coordinate systems. More...

Inheritance diagram for NURBS3DVolume:
[legend]
Collaboration diagram for NURBS3DVolume:
[legend]

Public Member Functions

 TypeName ("NURBS3DVolume")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, NURBS3DVolume, dictionary,(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors),(dict, mesh, computeParamCoors))
 
 NURBS3DVolume (const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
 Construct from dictionary. More...
 
 NURBS3DVolume (const NURBS3DVolume &)
 Construct as copy. More...
 
virtual ~NURBS3DVolume ()=default
 Destructor. More...
 
vector volumeDerivativeU (const scalar u, const scalar v, const scalar w) const
 Volume point derivative wrt u at point u,v,w. More...
 
vector volumeDerivativeV (const scalar u, const scalar v, const scalar w) const
 Volume point derivative wrt v at point u,v,w. More...
 
vector volumeDerivativeW (const scalar u, const scalar v, const scalar w) const
 Volume point derivative wrt w at point u,v,w. More...
 
tensor JacobianUVW (const vector &u) const
 Jacobian matrix wrt to the volume parametric coordinates. More...
 
scalar volumeDerivativeCP (const vector &u, const label cpI) const
 Volume point derivative wrt to control point cpI at point u,v,w. More...
 
vectorField computeControlPointSensitivities (const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
 
vectorField computeControlPointSensitivities (const volVectorField &faceSens, const labelList &sensitivityPatchIDs)
 
vectorField computeControlPointSensitivities (const boundaryVectorField &faceSens, const labelList &sensitivityPatchIDs)
 
vector computeControlPointSensitivities (const vectorField &faceSens, const label patchI, const label cpI)
 
tmp< tensorFielddndbBasedSensitivities (const label patchI, const label cpI, bool DimensionedNormalSens=true)
 
tmp< tensorFieldpatchDxDb (const label patchI, const label cpI)
 Get patch dx/db. More...
 
tmp< tensorFieldpatchDxDbFace (const label patchI, const label cpI)
 Get patch dx/db. More...
 
tmp< vectorFieldcoordinates (const vectorField &uVector) const
 
vector coordinates (const vector &uVector) const
 The same, for a specific point. More...
 
tmp< vectorFieldcomputeNewPoints (const vectorField &controlPointsMovement)
 Mesh movement based on given control point movement. More...
 
tmp< vectorFieldcomputeNewBoundaryPoints (const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
 Boundary mesh movement based on given control point movement. More...
 
void setControlPoints (const vectorField &newCps)
 Set new control points. More...
 
void boundControlPointMovement (vectorField &controlPointsMovement)
 
scalar computeMaxBoundaryDisplacement (const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
 Compute max. displacement at the boundary. More...
 
tmp< vectorFieldgetPointsInBox ()
 
const labelListgetMap ()
 Get map of points in box to mesh points. More...
 
const labelListgetReverseMap ()
 
const pointVectorFieldgetParametricCoordinates ()
 Get parametric coordinates. More...
 
tmp< pointTensorFieldgetDxDb (const label cpI)
 Get dxCartesiandb for a certain control point. More...
 
tmp< volTensorFieldgetDxCellsDb (const label cpI)
 Get dxCartesiandb for a certain control point on cells. More...
 
label nUSymmetry () const
 Get number of variables if CPs are moved symmetrically in U. More...
 
label nVSymmetry () const
 Get number of variables if CPs are moved symmetrically in V. More...
 
label nWSymmetry () const
 Get number of variables if CPs are moved symmetrically in W. More...
 
const wordname () const
 Get box name. More...
 
const boolListgetActiveCPs () const
 Which control points are active? More...
 
const boolListgetActiveDesignVariables () const
 Which design variables are active? More...
 
const vectorFieldgetControlPoints () const
 Get control points. More...
 
bool confineX1movement () const
 Get confine movements. More...
 
bool confineX2movement () const
 
bool confineX3movement () const
 
const NURBSbasisbasisU () const
 Get basis functions. More...
 
const NURBSbasisbasisV () const
 
const NURBSbasisbasisW () const
 
void writeCps (const string="cpsFile") const
 
void writeCpsInDict () const
 Write control points on the local coordinate system. More...
 
void write () const
 Write parametric coordinates. More...
 

Static Public Member Functions

static autoPtr< NURBS3DVolumeNew (const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
 Return a reference to the selected NURBS model. More...
 

Protected Types

typedef List< FixedList< bool, 3 > > boolListList3
 

Protected Member Functions

label getCPID (const label i, const label j, const label k) const
 Get control point ID from its I-J-K coordinates. More...
 
void findPointsInBox (const vectorField &meshPoints)
 Find points within control points box. More...
 
void computeParametricCoordinates (const vectorField &points)
 
void computeParametricCoordinates (tmp< vectorField > tPoints)
 
bool bound (vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999)
 Bound components to certain limits. More...
 
void determineActiveDesignVariablesAndPoints ()
 Create lists with active design variables and control points. More...
 
void confineBoundaryControlPoints ()
 Confine movement in boundary control points if necessary. More...
 
void continuityRealatedConfinement ()
 Confine control point movement to maintain user-defined continuity. More...
 
void confineControlPointsDirections ()
 Confine movement in all control points for user-defined directions. More...
 
void confineControlPoint (const label cpI)
 Confine all three movements for a prescribed control point. More...
 
void confineControlPoint (const label cpI, const FixedList< bool, 3 > &)
 Confine specific movements for a prescribed control point. More...
 
void makeFolders ()
 Create folders to store cps and derivatives. More...
 
virtual vector transformPointToCartesian (const vector &localCoordinates) const =0
 Transform a point from its coordinate system to a cartesian system. More...
 
virtual tensor transformationTensorDxDb (label globalPointIndex)=0
 
virtual void updateLocalCoordinateSystem (const vectorField &cartesianPoints)=0
 Update coordinates in the local system based on the cartesian points. More...
 

Protected Attributes

const fvMeshmesh_
 NURBS basis functions. More...
 
word name_
 
NURBSbasis basisU_
 
NURBSbasis basisV_
 
NURBSbasis basisW_
 
label maxIter_
 Max iterations for Newton-Raphson. More...
 
scalar tolerance_
 Tolerance for Newton-Raphson. More...
 
label nMaxBound_
 
vectorField cps_
 The volumetric B-Splines control points. More...
 
autoPtr< labelListmapPtr_
 Map of points-in-box to mesh points. More...
 
autoPtr< labelListreverseMapPtr_
 Map of mesh points to points-in-box. More...
 
autoPtr< pointVectorFieldparametricCoordinatesPtr_
 Parametric coordinates of pointsInBox. More...
 
vectorField localSystemCoordinates_
 Coordinates in the local system for which CPs are defined. More...
 
label confineX1movement_
 
label confineX2movement_
 
label confineX3movement_
 
label confineBoundaryControlPoints_
 
boolListList3 boundUMinCPs_
 Which movement components to freeze in each plane. More...
 
boolListList3 boundUMaxCPs_
 
boolListList3 boundVMinCPs_
 
boolListList3 boundVMaxCPs_
 
boolListList3 boundWMinCPs_
 
boolListList3 boundWMaxCPs_
 
boolList activeControlPoints_
 Which of the cps are moved in an optimisation. More...
 
boolList activeDesignVariables_
 Which design variables are changed in an optimisation. More...
 
string cpsFolder_
 Folder to store control points. More...
 
bool readStoredData_
 Read parametric coordinates from file if present in the folder. More...
 

Detailed Description

NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing support for different coordinate systems.

Reference:

    For a short introduction to a volumetric B-Splines morpher and its use
    in shape optimisation

        Papoutsis-Kiachagias, E., Magoulas, N., Mueller, J.,
        Othmer, C., & Giannakoglou, K. (2015).
        Noise reduction in car aerodynamics using a surrogate objective
        function and the continuous adjoint method with wall functions.
        Computers & Fluids, 122, 223-232.
        http://doi.org/10.1016/j.compfluid.2015.09.002
Source files

Definition at line 71 of file NURBS3DVolume.H.

Member Typedef Documentation

◆ boolListList3

typedef List<FixedList<bool, 3> > boolListList3
protected

Definition at line 77 of file NURBS3DVolume.H.

Constructor & Destructor Documentation

◆ NURBS3DVolume() [1/2]

NURBS3DVolume ( const dictionary dict,
const fvMesh mesh,
bool  computeParamCoors = true 
)

Construct from dictionary.

Definition at line 641 of file NURBS3DVolume.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, dictionary::get(), Foam::Info, IOobject::MUST_READ, IOobject::NO_WRITE, dictionary::readEntry(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ NURBS3DVolume() [2/2]

Construct as copy.

◆ ~NURBS3DVolume()

virtual ~NURBS3DVolume ( )
virtualdefault

Destructor.

Member Function Documentation

◆ getCPID()

Foam::label getCPID ( const label  i,
const label  j,
const label  k 
) const
protected

Get control point ID from its I-J-K coordinates.

Definition at line 50 of file NURBS3DVolume.C.

References k.

◆ findPointsInBox()

void findPointsInBox ( const vectorField meshPoints)
protected

Find points within control points box.

Definition at line 63 of file NURBS3DVolume.C.

References Field< Type >::component(), Foam::BitOps::count(), NURBS3DVolume::cps_, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::Info, NURBS3DVolume::mapPtr_, Foam::max(), Foam::min(), Foam::reduce(), NURBS3DVolume::reverseMapPtr_, List< T >::setSize(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ computeParametricCoordinates() [1/2]

void computeParametricCoordinates ( const vectorField points)
protected

Compute parametric coordinates in order to match a given set of coordinates, based on the cps of the class

Uses a Newton-Raphson loop. Argument is the points residing in the box

Definition at line 118 of file NURBS3DVolume.C.

References IOobject::AUTO_WRITE, Foam::bound(), Foam::cmptMag(), VectorSpace< Form, Cmpt, Ncmpts >::component(), coordinates(), Foam::diff(), Foam::dimless, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::Info, Foam::inv(), Foam::mag(), Foam::max(), Foam::min(), MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New(), IOobject::NO_READ, points, Foam::reduce(), tmp< T >::ref(), List< T >::setSize(), WarningInFunction, VectorSpace< Vector< scalar >, scalar, 3 >::zero, and Foam::Zero.

Here is the call graph for this function:

◆ computeParametricCoordinates() [2/2]

void computeParametricCoordinates ( tmp< vectorField tPoints)
protected

Definition at line 352 of file NURBS3DVolume.C.

References points.

◆ bound()

bool bound ( vector vec,
scalar  minValue = 1e-7,
scalar  maxValue = 0.999999 
)
protected

Bound components to certain limits.

Definition at line 362 of file NURBS3DVolume.C.

References maxValue, minValue, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ determineActiveDesignVariablesAndPoints()

void determineActiveDesignVariablesAndPoints ( )
protected

Create lists with active design variables and control points.

Definition at line 414 of file NURBS3DVolume.C.

References forAll.

◆ confineBoundaryControlPoints()

void confineBoundaryControlPoints ( )
protected

Confine movement in boundary control points if necessary.

Definition at line 446 of file NURBS3DVolume.C.

◆ continuityRealatedConfinement()

void continuityRealatedConfinement ( )
protected

Confine control point movement to maintain user-defined continuity.

Definition at line 492 of file NURBS3DVolume.C.

References forAll.

◆ confineControlPointsDirections()

void confineControlPointsDirections ( )
protected

Confine movement in all control points for user-defined directions.

Definition at line 587 of file NURBS3DVolume.C.

◆ confineControlPoint() [1/2]

void confineControlPoint ( const label  cpI)
protected

Confine all three movements for a prescribed control point.

Definition at line 598 of file NURBS3DVolume.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ confineControlPoint() [2/2]

void confineControlPoint ( const label  cpI,
const FixedList< bool, 3 > &  confineDirections 
)
protected

Confine specific movements for a prescribed control point.

Definition at line 617 of file NURBS3DVolume.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ makeFolders()

void makeFolders ( )
protected

Create folders to store cps and derivatives.

Definition at line 405 of file NURBS3DVolume.C.

References UPstream::master(), and Foam::mkDir().

Here is the call graph for this function:

◆ transformPointToCartesian()

virtual vector transformPointToCartesian ( const vector localCoordinates) const
protectedpure virtual

Transform a point from its coordinate system to a cartesian system.

Implemented in NURBS3DVolumeCylindrical, and NURBS3DVolumeCartesian.

◆ transformationTensorDxDb()

virtual tensor transformationTensorDxDb ( label  globalPointIndex)
protectedpure virtual

Transformation tensor for dxdb, from local coordinate system to cartesian

Implemented in NURBS3DVolumeCylindrical, and NURBS3DVolumeCartesian.

◆ updateLocalCoordinateSystem()

virtual void updateLocalCoordinateSystem ( const vectorField cartesianPoints)
protectedpure virtual

Update coordinates in the local system based on the cartesian points.

Implemented in NURBS3DVolumeCylindrical, and NURBS3DVolumeCartesian.

◆ TypeName()

TypeName ( "NURBS3DVolume"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
NURBS3DVolume  ,
dictionary  ,
(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors)  ,
(dict, mesh, computeParamCoors)   
)

◆ New()

Foam::autoPtr< Foam::NURBS3DVolume > New ( const dictionary dict,
const fvMesh mesh,
bool  computeParamCoors = true 
)
static

Return a reference to the selected NURBS model.

Definition at line 802 of file NURBS3DVolume.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::get(), Foam::Info, and mesh.

Referenced by volBSplinesBase::volBSplinesBase().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ volumeDerivativeU()

Foam::vector volumeDerivativeU ( const scalar  u,
const scalar  v,
const scalar  w 
) const

Volume point derivative wrt u at point u,v,w.

Definition at line 832 of file NURBS3DVolume.C.

References Foam::Zero.

◆ volumeDerivativeV()

Foam::vector volumeDerivativeV ( const scalar  u,
const scalar  v,
const scalar  w 
) const

Volume point derivative wrt v at point u,v,w.

Definition at line 868 of file NURBS3DVolume.C.

References Foam::Zero.

◆ volumeDerivativeW()

Foam::vector volumeDerivativeW ( const scalar  u,
const scalar  v,
const scalar  w 
) const

Volume point derivative wrt w at point u,v,w.

Definition at line 904 of file NURBS3DVolume.C.

References Foam::Zero.

◆ JacobianUVW()

Foam::tensor JacobianUVW ( const vector u) const

Jacobian matrix wrt to the volume parametric coordinates.

Definition at line 940 of file NURBS3DVolume.C.

References VectorSpace< Form, Cmpt, Ncmpts >::component(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), Vector< Cmpt >::z(), and Foam::Zero.

Here is the call graph for this function:

◆ volumeDerivativeCP()

Foam::scalar volumeDerivativeCP ( const vector u,
const label  cpI 
) const

Volume point derivative wrt to control point cpI at point u,v,w.

Scalar since in the local system!

Definition at line 969 of file NURBS3DVolume.C.

References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ computeControlPointSensitivities() [1/4]

Foam::vectorField computeControlPointSensitivities ( const pointVectorField pointSens,
const labelList sensitivityPatchIDs 
)

Control point sensitivities computed using point-based surface sensitivities

Definition at line 1002 of file NURBS3DVolume.C.

References forAll, Pstream::listCombineGather(), Pstream::listCombineScatter(), Foam::foamVersion::patch, and Foam::Zero.

Here is the call graph for this function:

◆ computeControlPointSensitivities() [2/4]

Foam::vectorField computeControlPointSensitivities ( const volVectorField faceSens,
const labelList sensitivityPatchIDs 
)

Control point sensitivities computed using face-based surface sensitivities

Definition at line 1053 of file NURBS3DVolume.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField().

Here is the call graph for this function:

◆ computeControlPointSensitivities() [3/4]

Foam::vectorField computeControlPointSensitivities ( const boundaryVectorField faceSens,
const labelList sensitivityPatchIDs 
)

Control point sensitivities computed using face-based surface sensitivities

Definition at line 1068 of file NURBS3DVolume.C.

References forAll, Pstream::listCombineGather(), Pstream::listCombineScatter(), deltaBoundary::makeFaceCentresAndAreas_d(), Foam::foamVersion::patch, face::points(), and Foam::Zero.

Here is the call graph for this function:

◆ computeControlPointSensitivities() [4/4]

Foam::vector computeControlPointSensitivities ( const vectorField faceSens,
const label  patchI,
const label  cpI 
)

Control point sensitivities computed using face-based surface sensitivities

Definition at line 1139 of file NURBS3DVolume.C.

References forAll, deltaBoundary::makeFaceCentresAndAreas_d(), Foam::foamVersion::patch, face::points(), Foam::reduce(), and Foam::Zero.

Here is the call graph for this function:

◆ dndbBasedSensitivities()

Foam::tmp< Foam::tensorField > dndbBasedSensitivities ( const label  patchI,
const label  cpI,
bool  DimensionedNormalSens = true 
)

Part of control point sensitivities related to the face normal variations

Definition at line 1199 of file NURBS3DVolume.C.

References forAll, deltaBoundary::makeFaceCentresAndAreas_d(), Foam::foamVersion::patch, face::points(), tmp< T >::ref(), polyPatch::start(), and Foam::Zero.

Here is the call graph for this function:

◆ patchDxDb()

Foam::tmp< Foam::tensorField > patchDxDb ( const label  patchI,
const label  cpI 
)

Get patch dx/db.

Definition at line 1268 of file NURBS3DVolume.C.

References forAll, tmp< T >::New(), Foam::foamVersion::patch, and Foam::Zero.

Here is the call graph for this function:

◆ patchDxDbFace()

Foam::tmp< Foam::tensorField > patchDxDbFace ( const label  patchI,
const label  cpI 
)

Get patch dx/db.

Definition at line 1307 of file NURBS3DVolume.C.

References forAll, deltaBoundary::makeFaceCentresAndAreas_d(), tmp< T >::New(), Foam::foamVersion::patch, face::points(), and Foam::Zero.

Here is the call graph for this function:

◆ coordinates() [1/2]

Foam::tmp< Foam::vectorField > coordinates ( const vectorField uVector) const

Compute cartesian coordinates based on control points and parametric coordinates

Definition at line 1365 of file NURBS3DVolume.C.

References forAll, tmp< T >::New(), nPoints, points, and Foam::Zero.

Here is the call graph for this function:

◆ coordinates() [2/2]

Foam::vector coordinates ( const vector uVector) const

The same, for a specific point.

Definition at line 1408 of file NURBS3DVolume.C.

References Vector< Cmpt >::x(), Vector< Cmpt >::y(), Vector< Cmpt >::z(), and Foam::Zero.

Here is the call graph for this function:

◆ computeNewPoints()

Foam::tmp< Foam::vectorField > computeNewPoints ( const vectorField controlPointsMovement)

Mesh movement based on given control point movement.

Definition at line 1445 of file NURBS3DVolume.C.

References coordinates(), DebugInfo, Foam::endl(), forAll, Foam::gMax(), Foam::mag(), and tmp< T >::ref().

Here is the call graph for this function:

◆ computeNewBoundaryPoints()

Foam::tmp< Foam::vectorField > computeNewBoundaryPoints ( const vectorField controlPointsMovement,
const labelList patchesToBeMoved 
)

Boundary mesh movement based on given control point movement.

Definition at line 1483 of file NURBS3DVolume.C.

References coordinates(), DebugInfo, Foam::endl(), Foam::gMax(), Foam::mag(), Foam::foamVersion::patch, and tmp< T >::ref().

Here is the call graph for this function:

◆ setControlPoints()

void setControlPoints ( const vectorField newCps)

Set new control points.

New values should be on the coordinates system original CPs were defined

Definition at line 1536 of file NURBS3DVolume.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ boundControlPointMovement()

void boundControlPointMovement ( vectorField controlPointsMovement)

Bound control points movement in the boundary control points and in certain directions if needed

Definition at line 1550 of file NURBS3DVolume.C.

References forAll, and Foam::Zero.

◆ computeMaxBoundaryDisplacement()

Foam::scalar computeMaxBoundaryDisplacement ( const vectorField controlPointsMovement,
const labelList patchesToBeMoved 
)

Compute max. displacement at the boundary.

Definition at line 1573 of file NURBS3DVolume.C.

References coordinates(), Foam::mag(), Foam::max(), Foam::foamVersion::patch, Foam::reduce(), and Foam::Zero.

Here is the call graph for this function:

◆ getPointsInBox()

Foam::tmp< Foam::vectorField > getPointsInBox ( )

Get points of mesh which recide within the control points box

Definition at line 1622 of file NURBS3DVolume.C.

◆ getMap()

const Foam::labelList & getMap ( )

Get map of points in box to mesh points.

Definition at line 1637 of file NURBS3DVolume.C.

◆ getReverseMap()

const Foam::labelList & getReverseMap ( )

Get map of mesh points to points in box. Return -1 if point is outside the box

Definition at line 1648 of file NURBS3DVolume.C.

◆ getParametricCoordinates()

const Foam::pointVectorField & getParametricCoordinates ( )

Get parametric coordinates.

Definition at line 1659 of file NURBS3DVolume.C.

◆ getDxDb()

Foam::tmp< Foam::pointTensorField > getDxDb ( const label  cpI)

Get dxCartesiandb for a certain control point.

Definition at line 1677 of file NURBS3DVolume.C.

References Foam::dimless, MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New(), IOobject::NO_READ, IOobject::NO_WRITE, tmp< T >::ref(), and Foam::Zero.

Here is the call graph for this function:

◆ getDxCellsDb()

Foam::tmp< Foam::volTensorField > getDxCellsDb ( const label  cpI)

Get dxCartesiandb for a certain control point on cells.

Definition at line 1721 of file NURBS3DVolume.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), deltaBoundary::cellCenters_d(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimless, forAll, IOobject::NO_READ, IOobject::NO_WRITE, Foam::foamVersion::patch, tmp< T >::ref(), and Foam::Zero.

Here is the call graph for this function:

◆ nUSymmetry()

Foam::label nUSymmetry ( ) const

Get number of variables if CPs are moved symmetrically in U.

Definition at line 1789 of file NURBS3DVolume.C.

◆ nVSymmetry()

Foam::label nVSymmetry ( ) const

Get number of variables if CPs are moved symmetrically in V.

Definition at line 1804 of file NURBS3DVolume.C.

◆ nWSymmetry()

Foam::label nWSymmetry ( ) const

Get number of variables if CPs are moved symmetrically in W.

Definition at line 1819 of file NURBS3DVolume.C.

◆ name()

const Foam::word & name ( ) const
inline

Get box name.

Definition at line 34 of file NURBS3DVolumeI.H.

References NURBS3DVolume::name_.

◆ getActiveCPs()

const Foam::boolList & getActiveCPs ( ) const
inline

Which control points are active?

A control point is active if at least one component can move

Definition at line 40 of file NURBS3DVolumeI.H.

◆ getActiveDesignVariables()

const Foam::boolList & getActiveDesignVariables ( ) const
inline

Which design variables are active?

Numbering is (X1,Y1,Z1), (X2,Y2,Z2) ...

Definition at line 47 of file NURBS3DVolumeI.H.

◆ getControlPoints()

const Foam::vectorField & getControlPoints ( ) const
inline

Get control points.

Definition at line 53 of file NURBS3DVolumeI.H.

◆ confineX1movement()

bool confineX1movement ( ) const
inline

Get confine movements.

Definition at line 59 of file NURBS3DVolumeI.H.

◆ confineX2movement()

bool confineX2movement ( ) const
inline

Definition at line 65 of file NURBS3DVolumeI.H.

◆ confineX3movement()

bool confineX3movement ( ) const
inline

Definition at line 71 of file NURBS3DVolumeI.H.

◆ basisU()

const Foam::NURBSbasis & basisU ( ) const
inline

Get basis functions.

Definition at line 77 of file NURBS3DVolumeI.H.

◆ basisV()

const Foam::NURBSbasis & basisV ( ) const
inline

Definition at line 83 of file NURBS3DVolumeI.H.

◆ basisW()

const Foam::NURBSbasis & basisW ( ) const
inline

Definition at line 89 of file NURBS3DVolumeI.H.

◆ writeCps()

void writeCps ( const string  fileName = "cpsFile") const

Write control points on a cartesian coordinates system for visualization

Definition at line 1834 of file NURBS3DVolume.C.

References Foam::endl(), forAll, Foam::Info, UPstream::master(), and Foam::Zero.

Here is the call graph for this function:

◆ writeCpsInDict()

void writeCpsInDict ( ) const

Write control points on the local coordinate system.

For continuation

Definition at line 1876 of file NURBS3DVolume.C.

References dictionary::add(), IOstreamOption::ASCII, IOstreamOption::currentVersion, UPstream::master(), IOobject::NO_READ, and IOobject::NO_WRITE.

Here is the call graph for this function:

◆ write()

void write ( ) const

Write parametric coordinates.

Definition at line 1911 of file NURBS3DVolume.C.

Member Data Documentation

◆ mesh_

const fvMesh& mesh_
protected

NURBS basis functions.

Definition at line 80 of file NURBS3DVolume.H.

◆ name_

word name_
protected

Definition at line 81 of file NURBS3DVolume.H.

Referenced by NURBS3DVolume::name().

◆ basisU_

NURBSbasis basisU_
protected

Definition at line 82 of file NURBS3DVolume.H.

◆ basisV_

NURBSbasis basisV_
protected

Definition at line 83 of file NURBS3DVolume.H.

◆ basisW_

NURBSbasis basisW_
protected

Definition at line 84 of file NURBS3DVolume.H.

◆ maxIter_

label maxIter_
protected

Max iterations for Newton-Raphson.

Definition at line 87 of file NURBS3DVolume.H.

◆ tolerance_

scalar tolerance_
protected

Tolerance for Newton-Raphson.

Definition at line 90 of file NURBS3DVolume.H.

◆ nMaxBound_

label nMaxBound_
protected

How many times to bound parametric coordinates until deciding it is outside the box

Definition at line 94 of file NURBS3DVolume.H.

◆ cps_

vectorField cps_
protected

The volumetric B-Splines control points.

Definition at line 97 of file NURBS3DVolume.H.

Referenced by NURBS3DVolume::findPointsInBox().

◆ mapPtr_

autoPtr<labelList> mapPtr_
protected

Map of points-in-box to mesh points.

Definition at line 100 of file NURBS3DVolume.H.

Referenced by NURBS3DVolume::findPointsInBox().

◆ reverseMapPtr_

autoPtr<labelList> reverseMapPtr_
protected

Map of mesh points to points-in-box.

Return -1 if point not whithin the box

Definition at line 104 of file NURBS3DVolume.H.

Referenced by NURBS3DVolume::findPointsInBox().

◆ parametricCoordinatesPtr_

autoPtr<pointVectorField> parametricCoordinatesPtr_
protected

Parametric coordinates of pointsInBox.

Definition at line 107 of file NURBS3DVolume.H.

◆ localSystemCoordinates_

vectorField localSystemCoordinates_
protected

Coordinates in the local system for which CPs are defined.

Definition at line 110 of file NURBS3DVolume.H.

◆ confineX1movement_

label confineX1movement_
protected

Confine movement in certain directions and control points. Refers to the local system

Definition at line 114 of file NURBS3DVolume.H.

◆ confineX2movement_

label confineX2movement_
protected

Definition at line 116 of file NURBS3DVolume.H.

◆ confineX3movement_

label confineX3movement_
protected

Definition at line 118 of file NURBS3DVolume.H.

◆ confineBoundaryControlPoints_

label confineBoundaryControlPoints_
protected

Definition at line 120 of file NURBS3DVolume.H.

◆ boundUMinCPs_

boolListList3 boundUMinCPs_
protected

Which movement components to freeze in each plane.

Definition at line 123 of file NURBS3DVolume.H.

◆ boundUMaxCPs_

boolListList3 boundUMaxCPs_
protected

Definition at line 125 of file NURBS3DVolume.H.

◆ boundVMinCPs_

boolListList3 boundVMinCPs_
protected

Definition at line 127 of file NURBS3DVolume.H.

◆ boundVMaxCPs_

boolListList3 boundVMaxCPs_
protected

Definition at line 129 of file NURBS3DVolume.H.

◆ boundWMinCPs_

boolListList3 boundWMinCPs_
protected

Definition at line 131 of file NURBS3DVolume.H.

◆ boundWMaxCPs_

boolListList3 boundWMaxCPs_
protected

Definition at line 133 of file NURBS3DVolume.H.

◆ activeControlPoints_

boolList activeControlPoints_
protected

Which of the cps are moved in an optimisation.

Definition at line 136 of file NURBS3DVolume.H.

◆ activeDesignVariables_

boolList activeDesignVariables_
protected

Which design variables are changed in an optimisation.

Definition at line 139 of file NURBS3DVolume.H.

◆ cpsFolder_

string cpsFolder_
protected

Folder to store control points.

Definition at line 142 of file NURBS3DVolume.H.

◆ readStoredData_

bool readStoredData_
protected

Read parametric coordinates from file if present in the folder.

Definition at line 145 of file NURBS3DVolume.H.


The documentation for this class was generated from the following files: