59#ifndef topODesignVariables_H
60#define topODesignVariables_H
75class topODesignVariables
77 public topOVariablesBase,
78 public designVariables
108 const label fluidID = 0
117 const label fluidID = 0,
118 const bool activeIO =
false
126 const label fluidID = 0,
127 const bool setIOValues =
false
139 topODesignVariables(
const topODesignVariables&) =
delete;
142 void operator=(
const topODesignVariables&) =
delete;
191 const word& interpolationFieldName =
"beta"
202 const word& interpolationFieldName =
"beta"
214 const word& designVariablesName,
215 const word& interpolationFieldName =
"beta"
230 const word& designVariablesName
236 const word& interpolationFieldName,
243 const word& interpolationFieldName,
A field of fields is a PtrList of fields with reference counting.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
SubField is a Field obtained as a section of another Field, without its own allocation....
Abstract base class for adjoint-based sensitivities.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
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.
Base class for primal solvers.
A class for managing temporary objects.
TypeName("topO")
Runtime type information.
virtual void updateField(const scalarField &correction, const label fluidID=0)
Update the design variables given their correction.
virtual const scalarField & interpolationField(const word &interpolationFieldName="beta") const
Return interpolant.
fieldRegularisation & getRegularisation()
Get access to the regularisation object.
void readField(const word &name, const label fluidID=0, const bool setIOValues=false)
Read field with (path of) the design variables and store input in the design variables list with opti...
fieldRegularisation regularisation_
Mechanism to regularise the field of design variables.
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
virtual bool writeData(Ostream &) const
The writeData function required by the regIOobject write operation.
virtual void initialize()
Part of the constructor initialisation.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &sens)
Assemble sensitivity derivatives, by combining the part related to the primal and adjoint solution wi...
virtual tmp< scalarField > penaltySensitivities(const word &interpolationFieldName, const topOInterpolationFunction &interpolationFunc) const
Return the penalty term derivative.
bool addFvOptions_
Add the fvOptions necessary for topO automatically.
virtual void nullifyInSolidSensitivities(scalarField &sens, const topOInterpolationFunction &interpolationFunc, const word &designVariablesName) const
Nullify given field in the solid domain.
virtual void addFvOptions(const PtrList< primalSolver > &primalSolver, const PtrList< adjointSolverManager > &adjointSolverManagers)
Automatically add fvOptions depending on the design variables to the primal and adjoint solvers.
SubField< scalar > alpha_
A subfield of the design variables corresponding to the porosity field.
static autoPtr< topODesignVariables > New(fvMesh &mesh, const dictionary &dict)
Construct and return the selected design variables.
virtual void setActiveDesignVariables(const label fluidID=0, const bool activeIO=false)
Set active design variables.
bool writeAllFields_
Write all fields related to imposition of the Brinkman penalisation (i.e. design variables,...
virtual void interpolate(volScalarField &field, const topOInterpolationFunction &interpolationFunc, const FieldField< Field, scalar > &fluidValues, const scalarField &solidValues, const label fieldi, const word &interpolationFieldName="beta") const
Interpolate between the given field and solid values.
virtual void interpolationSensitivities(scalarField &sens, const topOInterpolationFunction &interpolationFunc, const FieldField< Field, scalar > &fluidValues, const scalarField &solidValues, const label fieldi, const word &designVariablesName, const word &interpolationFieldName="beta") const
Post-processing sensitivities due to interpolations based on the indicator fields.
virtual const volScalarField & beta() const
Get the indicator field.
virtual void setInitialValues()
Set initial values of the design variables.
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
virtual ~topODesignVariables()=default
Destructor.
void applyFixedValues()
Apply fixed values in certain zones.
virtual tmp< scalarField > penalty(const word &interpolationFieldName, const topOInterpolationFunction &interpolationFunc) const
Return the Brinkman penalisation term.
virtual void nullifyInSolid(scalarField &field, const topOInterpolationFunction &interpolationFunc) const
Nullify given field in the solid domain.
virtual void writeDesignVars()
Write porosity field to file.
A class for handling words, derived from Foam::string.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
PtrList< adjointSolverManager > & adjointSolverManagers