45#ifndef volBSplinesBase_H
46#define volBSplinesBase_H
62 public MeshObject<fvMesh, UpdateableMeshObject, volBSplinesBase>
77 volBSplinesBase(
const volBSplinesBase&) =
delete;
80 void operator=(
const volBSplinesBase&) =
delete;
MeshObject(const fvMesh &mesh)
const fvMesh & mesh() const noexcept
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Mesh data needed to do the Finite Volume discretisation.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
virtual bool movePoints()
Dummy function required by MeshObject.
virtual ~volBSplinesBase()=default
Destructor.
void moveControlPoints(const vectorField &controlPointsMovement)
Move control points. No effect on mesh.
const vectorField & getControlPoints(const label &iNURB) const
Get reference to control points.
const PtrList< NURBS3DVolume > & boxes() const
Get const reference to the vol. B-splines boxes.
vectorField getAllControlPoints() const
Get control points from all boxes.
const labelList & getActiveDesignVariables() const
Get active design variables.
void writeControlPoints() const
Write control points to constant and optimisation folders.
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Get max boundary displacement for a given control-points movement.
labelList getStartCpID() const
Get start CP ID for each box.
labelList getStartVarID() const
Get start CP ID for each box.
label findBoxID(const label cpI) const
Find box of certain control point.
NURBS3DVolume & boxRef(const label boxI)
Get non-const reference to a specific box.
PtrList< NURBS3DVolume > & boxesRef()
Get non-const reference to the vol. B-splines boxes.
TypeName("volBSplinesBase")
Runtime type information.
labelList activeDesignVariables_
Active design variables numbering for all boxes.
const NURBS3DVolume & box(const label boxI) const
Get const reference to a specific box.
label getTotalControlPointsNumber() const
Get cumulative number of control points from all boxes.
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement.
PtrList< NURBS3DVolume > volume_
List with volumetric B-splines boxes.
virtual void updateMesh(const mapPolyMesh &)
Dummy function required by MeshObject.
tmp< vectorField > computeBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Get the updated boundary points only.
label getNumberOfBoxes() const
Get number of boxes.
Vector< label > decomposeDV(const label dvI) const
From design variable ID, return boxID, cpID and direction.
List< label > labelList
A List of labels.
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.