38#ifndef volumetricBSplinesDesignVariables_H
39#define volumetricBSplinesDesignVariables_H
54class volumetricBSplinesDesignVariables
56 public shapeDesignVariables,
112 const word& boundsName,
130 volumetricBSplinesDesignVariables
132 const volumetricBSplinesDesignVariables&
136 void operator=(
const volumetricBSplinesDesignVariables&) =
delete;
148 volumetricBSplinesDesignVariables
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Abstract base class for adjoint-based sensitivities.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
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.
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
virtual void evolveNumber()
For design variables with a dynamic character (i.e. changing.
virtual tmp< vectorField > dSdb(const label patchI, const label varID) const
Get dSdb for given design variable and patch.
virtual ~volumetricBSplinesDesignVariables()=default
Destructor.
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
virtual void resetDesignVariables()
Reset to starting point of line search.
const volBSplinesBase & getVolBSplinesBase() const
Constant access to the volBSplinesBase object.
bool updateBounds() const
Are bounds to be updated after each optmisation cycle.
void controlPointsToDesignVariables()
Set the design variables based on the current control points.
void setDisplacement(const vectorField &cpMovement)
Set the field driving the displacement method.
virtual bool writeData(Ostream &os) const
Write fields to support continuation.
void readBounds(autoPtr< scalar > lowerBoundPtr=nullptr, autoPtr< scalar > upperBoundPtr=nullptr)
Read bounds for design variables, if present.
void designVariablesToControlPoints()
Set control points based on current design variables values.
void setActiveDesignVariables()
Set IDs of active design variables.
TypeName("volumetricBSplines")
Runtime type information.
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
bool nonOverlappingCPs() const
Are control points non-overlapping.
void writeBounds(const scalarField &bounds, const word &name) const
Write current bounds to file.
bool updateBounds_
Should the bounds of the non-overlapping control points be updated in each optimisation cycle.
tmp< vectorField > controlPointMovement(const scalarField &correction)
Transform correction of design variables to control points movement.
volBSplinesBase & getVolBSplinesBase()
Non-constant access to the volBSplinesBase object.
virtual tmp< vectorField > dxdbVol(const label varID) const
Get dxdb for all mesh points.
virtual const labelList & activeSensitivities() const
Components of the active control points.
bool nonOverlappingCPs_
Should the control points be non-overlapping.
volBSplinesBase & volBSplinesBase_
Reference to underlaying volumetric B-Splines morpher.
virtual label sensSize() const
Size of the active control points.
autoPtr< morphingBoxConstraint > constraint_
Constraint imposed on the movement of the design variables.
virtual tmp< vectorField > dxdbFace(const label patchI, const label varID) const
Get dxdb for given design variable and patch.
virtual tmp< vectorField > dndb(const label patchI, const label varID) const
Get dndb for given design variable and patch.
virtual tmp< volVectorField > dCdb(const label varID) const
Get dCdb for given design variable.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Assemble the sensitivity derivatives, by also applying possible constraints.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
Namespace for bounding specifications. At the moment, mostly for tables.
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.