53#ifndef fieldRegularisation_H
54#define fieldRegularisation_H
73class fieldRegularisation
127 fieldRegularisation(
const fieldRegularisation&);
130 void operator=(
const fieldRegularisation&);
184 const bool isTopoField,
193 const bool isTopoField
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,...
bool growFromWalls() const
Grow beta from walls.
bool growFromWalls_
Whether to apply a fixedValue BC or zeroGradient one to alphaTilda, when regularisation is performed.
autoPtr< regularisationRadius > radius_
Smoothing radius.
bool regularise_
Perform regulaisation on alpha before inputing it on beta?
virtual void updateBeta()
Update the beta field.
TypeName("fieldRegularisation")
Runtime type information.
bool project_
Perform the projection (sharpening) step?
volScalarField beta_
Beta is the field used for all interpolations between fluid and solid in topology optimisation.
autoPtr< topOInterpolationFunction > sharpenFunction_
Function used to sharpen the field after regularisation.
void postProcessSens(scalarField &sens)
Update part of fieldRegularisation to the sensitivitiy derivatives.
virtual ~fieldRegularisation()=default
Destructor.
autoPtr< volScalarField > alphaTilda_
The regularised alpha field, if regulatisation is performed.
void regularise(const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius)
Regularise an externally provided radius.
const volScalarField & beta() const
Return beta field.
bool shouldRegularise() const
Should regularisation be executed.
autoPtr< regularisationPDE > regularisationPDE_
PDE used for the regularisation.
const topOZones & zones_
Cell zones related to topology optimisation.
const scalarField & betaArg_
Argument of the beta field.
const scalarField & alpha_
Alpha field (design variables of topology optimisation).
Mesh data needed to do the Finite Volume discretisation.
Base class for selecting the regulatisation radius.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.