54#ifndef NURBS3DVolume_H
55#define NURBS3DVolume_H
202 const vector& localCoordinates
232 bool computeParamCoors
244 bool computeParamCoors =
true
258 bool computeParamCoors =
true
350 bool DimensionedNormalSens =
true
388 const bool moveCPs =
true
395 label
getCPID(
const label i,
const label j,
const label
k)
const;
398 void getIJK(label& i, label& j, label&
k,
const label cpID)
const;
Useful typenames for fields defined only at the boundaries.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
bool readStoredData_
Read parametric coordinates from file if present in the folder.
void getIJK(label &i, label &j, label &k, const label cpID) const
Get I-J-K coordinates from control point ID.
label maxIter_
Max iterations for Newton-Raphson.
virtual void updateLocalCoordinateSystem(const vectorField &cartesianPoints)=0
Update coordinates in the local system based on the cartesian points.
boolVectorList confineWMinCPs_
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
boolVectorList confineVMinCPs_
tmp< tensorField > patchDxDbFace(const label patchI, const label cpI)
Get patch dx/db.
void computeParametricCoordinates(const vectorField &points)
Compute parametric coordinates in order to match a given set of coordinates, based on the cps of the ...
label nUSymmetry() const
Get number of variables if CPs are moved symmetrically in U.
void confineInertControlPoints()
Deactivate control points if they affect no mesh point.
autoPtr< labelList > mapPtr_
Map of points-in-box to mesh points.
declareRunTimeSelectionTable(autoPtr, NURBS3DVolume, dictionary,(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors),(dict, mesh, computeParamCoors))
vector volumeDerivativeU(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt u at point u,v,w.
void confineBoundaryControlPoints()
Confine movement in boundary control points if necessary.
const labelList & getMap()
Get map of points in box to mesh points.
scalar tolerance_
Tolerance for Newton-Raphson.
virtual bool writeData(Ostream &os) const
Write the control points to support restart.
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Compute max. displacement at the boundary.
boolVectorList confineVMaxCPs_
void findPointsInBox(const vectorField &meshPoints)
Find points within control points box.
tmp< vectorField > computeNewBoundaryPoints(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved, const bool moveCPs=true)
Boundary mesh movement based on given control point movement.
boolVectorList confineWMaxCPs_
const boolList & getActiveDesignVariables() const
Which design variables are active?
vectorField localSystemCoordinates_
Coordinates in the local system for which CPs are defined.
boolVectorList confineUMaxCPs_
vectorField computeControlPointSensitivities(const pointVectorField &pointSens, const labelList &sensitivityPatchIDs)
Control point sensitivities computed using point-based surface sensitivities.
boolVectorList confineUMinCPs_
Which movement components to freeze in each plane.
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
tmp< volTensorField > getDxCellsDb(const label cpI)
Get dxCartesiandb for a certain control point on cells.
virtual tensor transformationTensorDxDb(label globalPointIndex)=0
Transformation tensor for dxdb, from local coordinate system to cartesian.
autoPtr< labelList > reverseMapPtr_
Map of mesh points to points-in-box.
label nVSymmetry() const
Get number of variables if CPs are moved symmetrically in V.
tmp< pointTensorField > getDxDb(const label cpI)
Get dxCartesiandb for a certain control point.
NURBS3DVolume(const NURBS3DVolume &)
Construct as copy.
virtual ~NURBS3DVolume()=default
Destructor.
bool confineVMovement() const
bool bound(vector &vec, scalar minValue=1e-7, scalar maxValue=0.999999) const
Bound components to certain limits.
void makeFolders()
Create folders to store cps and derivatives.
void writeCps(const fileName &="cpsFile", const bool transform=true) const
Write control points on a cartesian coordinates system for visualization.
NURBS3DVolume(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Construct from dictionary.
boolList activeDesignVariables_
Which design variables are changed in an optimisation.
tmp< tensorField > patchDxDb(const label patchI, const label cpI)
Get patch dx/db.
tmp< vectorField > getPointsInBox()
Get mesh points that reside within the control points box.
vector volumeDerivativeV(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt v at point u,v,w.
virtual vector transformPointToCartesian(const vector &localCoordinates) const =0
Transform a point from its coordinate system to a cartesian system.
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool DimensionedNormalSens=true)
Part of control point sensitivities related to the face normal variations.
string cpsFolder_
Folder to store control points.
NURBSbasis basisU_
NURBS basis functions.
label confineUMovement_
Confine movement in certain directions and control points. Refers to the local system.
void confineControlPoint(const label cpI)
Confine all three movements for a prescribed control point.
bool confineWMovement() const
const dictionary & dict() const
Get dictionary.
Vector< label > nSymmetry() const
Get number of variables per direction, if CPs are moved symmetrically.
vectorField cps_
The volumetric B-Splines control points.
const NURBSbasis & basisV() const
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement in the boundary control points and in certain directions if needed.
void continuityRealatedConfinement()
Confine control point movement to maintain user-defined continuity.
label getCPID(const label i, const label j, const label k) const
Get control point ID from its I-J-K coordinates.
autoPtr< pointVectorField > parametricCoordinatesPtr_
Parametric coordinates of pointsInBox.
scalar volumeDerivativeCP(const vector &u, const label cpI) const
Volume point derivative wrt to control point cpI at point u,v,w.
label confineBoundaryControlPoints_
label nMaxBound_
How many times to bound parametric coordinates until deciding it is outside the box.
const NURBSbasis & basisW() const
const fvMesh & mesh() const
Get mesh.
label nWSymmetry() const
Get number of variables if CPs are moved symmetrically in W.
void writeParamCoordinates() const
Write parametric coordinates.
void determineActiveDesignVariablesAndPoints()
Create lists with active design variables and control points.
Vector< label > nCPsPerDirection() const
Get number of control points per direction.
tmp< vectorField > coordinates(const vectorField &uVector) const
Compute cartesian coordinates based on control points and parametric coordinates.
tmp< vectorField > computeNewPoints(const vectorField &controlPointsMovement)
Mesh movement based on given control point movement.
const word & name() const
Get box name.
List< boolVector > boolVectorList
bool confineUMovement() const
Get confine movements.
const boolList & getActiveCPs() const
Which control points are active?
tensor JacobianUVW(const vector &u) const
Jacobian matrix wrt to the volume parametric coordinates.
void setControlPoints(const vectorField &newCps)
Set new control points.
boolList activeControlPoints_
Which of the cps are moved in an optimisation.
const labelList & getReverseMap()
Get map of mesh points to points in box. Return -1 if point is outside the box.
TypeName("NURBS3DVolume")
Runtime type information.
const NURBSbasis & basisU() const
Get basis functions.
const vectorField & getControlPoints() const
Get control points.
vector volumeDerivativeW(const scalar u, const scalar v, const scalar w) const
Volume point derivative wrt w at point u,v,w.
void confineControlPointsDirections()
Confine movement in all control points for user-defined directions.
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Specialized bundling of boolean values as a vector of 3 components, element access using x(),...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A class for handling file names.
Mesh data needed to do the Finite Volume discretisation.
localIOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
volVectorField::Boundary boundaryVectorField
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.