Go to the documentation of this file.
54 #ifndef NURBS3DVolume_H
55 #define NURBS3DVolume_H
196 const vector& localCoordinates
226 bool computeParamCoors
238 bool computeParamCoors =
true
252 bool computeParamCoors =
true
337 bool DimensionedNormalSens =
true
381 label
getCPID(
const label i,
const label j,
const label
k)
const;
void computeParametricCoordinates(const vectorField &points)
const dictionary & dict() const
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
A class for handling words, derived from Foam::string.
const NURBSbasis & basisV() const
bool confineWMovement() const
A class for handling file names.
const labelList & getMap()
vector volumeDerivativeU(const scalar u, const scalar v, const scalar w) const
boolList activeDesignVariables_
A class for managing temporary objects.
boolList activeControlPoints_
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
boolVectorList confineVMaxCPs_
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
vector volumeDerivativeW(const scalar u, const scalar v, const scalar w) const
autoPtr< pointVectorField > parametricCoordinatesPtr_
virtual ~NURBS3DVolume()=default
vectorField localSystemCoordinates_
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
bool bound(vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999)
void confineControlPoint(const label cpI)
dimensionSet transform(const dimensionSet &ds)
virtual void updateLocalCoordinateSystem(const vectorField &cartesianPoints)=0
void confineControlPointsDirections()
Specialized bundling of boolean values as a vector of 3 components, element access using x(),...
void continuityRealatedConfinement()
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
autoPtr< labelList > reverseMapPtr_
label getCPID(const label i, const label j, const label k) const
const pointVectorField & getParametricCoordinates()
void confineBoundaryControlPoints()
TypeName("NURBS3DVolume")
List< boolVector > boolVectorList
Generic templated field type.
bool confineVMovement() const
void findPointsInBox(const vectorField &meshPoints)
const boolList & getActiveDesignVariables() const
boolVectorList confineWMaxCPs_
void determineActiveDesignVariablesAndPoints()
autoPtr< labelList > mapPtr_
scalar volumeDerivativeCP(const vector &u, const label cpI) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
void writeCpsInDict() const
Mesh data needed to do the Finite Volume discretisation.
label confineBoundaryControlPoints_
tmp< vectorField > getPointsInBox()
const NURBSbasis & basisW() const
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
vector volumeDerivativeV(const scalar u, const scalar v, const scalar w) const
boolVectorList confineVMinCPs_
virtual tensor transformationTensorDxDb(label globalPointIndex)=0
const vectorField & getControlPoints() const
void setControlPoints(const vectorField &newCps)
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
Useful typenames for fields defined only at the boundaries.
const boolList & getActiveCPs() const
boolVectorList confineUMinCPs_
virtual vector transformPointToCartesian(const vector &localCoordinates) const =0
const labelList & getReverseMap()
void boundControlPointMovement(vectorField &controlPointsMovement)
tensor JacobianUVW(const vector &u) const
const dimensionedScalar e
tmp< vectorField > coordinates(const vectorField &uVector) const
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
const word & name() const
bool confineUMovement() const
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
void writeCps(const fileName &="cpsFile", const bool transform=true) const
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
tmp< pointTensorField > getDxDb(const label cpI)
const fvMesh & mesh() const
tmp< volTensorField > getDxCellsDb(const label cpI)
const NURBSbasis & basisU() const
Generic GeometricField class.
tmp< tensorField > patchDxDb(const label patchI, const label cpI)
tmp< vectorField > computeNewPoints(const vectorField &controlPointsMovement)
boolVectorList confineWMinCPs_
boolVectorList confineUMaxCPs_
declareRunTimeSelectionTable(autoPtr, NURBS3DVolume, dictionary,(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors),(dict, mesh, computeParamCoors))