#include <fieldRegularisation.H>

Public Member Functions | |
| TypeName ("fieldRegularisation") | |
| Runtime type information. | |
| fieldRegularisation (fvMesh &mesh, const scalarField &alpha, const topOZones &zones, const dictionary &dict) | |
| Construct from components. | |
| virtual | ~fieldRegularisation ()=default |
| Destructor. | |
| const volScalarField & | beta () const |
| Return beta field. | |
| bool | growFromWalls () const |
| Grow beta from walls. | |
| bool | shouldRegularise () const |
| Should regularisation be executed. | |
| virtual void | updateBeta () |
| Update the beta field. | |
| void | regularise (const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius) |
| Regularise an externally provided radius. | |
| void | regularise (const scalarField &source, scalarField &result, const bool isTopoField) |
| Regularise field (or sensitivities) using a Laplace PDE. | |
| void | postProcessSens (scalarField &sens) |
| Update part of fieldRegularisation to the sensitivitiy derivatives. | |
Protected Attributes | |
| fvMesh & | mesh_ |
| dictionary | dict_ |
| const topOZones & | zones_ |
| Cell zones related to topology optimisation. | |
| bool | regularise_ |
| Perform regulaisation on alpha before inputing it on beta? | |
| bool | project_ |
| Perform the projection (sharpening) step? | |
| autoPtr< regularisationRadius > | radius_ |
| Smoothing radius. | |
| const scalarField & | alpha_ |
| Alpha field (design variables of topology optimisation). | |
| autoPtr< volScalarField > | alphaTilda_ |
| The regularised alpha field, if regulatisation is performed. | |
| autoPtr< topOInterpolationFunction > | sharpenFunction_ |
| Function used to sharpen the field after regularisation. | |
| autoPtr< regularisationPDE > | regularisationPDE_ |
| PDE used for the regularisation. | |
| const scalarField & | betaArg_ |
| Argument of the beta field. | |
| bool | growFromWalls_ |
| Whether to apply a fixedValue BC or zeroGradient one to alphaTilda, when regularisation is performed. | |
| volScalarField | beta_ |
| Beta is the field used for all interpolations between fluid and solid in topology optimisation. | |
Definition at line 68 of file fieldRegularisation.H.
| fieldRegularisation | ( | fvMesh & | mesh, |
| const scalarField & | alpha, | ||
| const topOZones & | zones, | ||
| const dictionary & | dict ) |
Construct from components.
Definition at line 36 of file fieldRegularisation.C.
References alpha, alpha_, alphaTilda_, beta_, betaArg_, DebugInfo, dict, dict_, Foam::dimless, Foam::dimTime, Foam::endl(), Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, growFromWalls_, mesh, mesh_, Foam::New(), Foam::nl, project_, radius_, regularisationPDE_, regularise_, sharpenFunction_, timeName, Foam::Zero, and zones_.

|
virtualdefault |
Destructor.
|
inline |
|
inline |
Grow beta from walls.
Definition at line 197 of file fieldRegularisation.H.
References growFromWalls_.
|
inline |
Should regularisation be executed.
Definition at line 205 of file fieldRegularisation.H.
References regularise_.
|
virtual |
Update the beta field.
Performs regulaisation of alpha_, if necessary
Definition at line 102 of file fieldRegularisation.C.
References alpha_, alphaTilda_, beta_, betaArg_, project_, regularise(), regularise_, and sharpenFunction_.

| void regularise | ( | const scalarField & | source, |
| scalarField & | result, | ||
| const bool | isTopoField, | ||
| const regularisationRadius & | radius ) |
Regularise an externally provided radius.
Definition at line 122 of file fieldRegularisation.C.
References alphaTilda_, regularisationPDE_, and regularise().
Referenced by postProcessSens(), regularise(), regularise(), and updateBeta().


| void regularise | ( | const scalarField & | source, |
| scalarField & | result, | ||
| const bool | isTopoField ) |
Regularise field (or sensitivities) using a Laplace PDE.
Definition at line 135 of file fieldRegularisation.C.
References alphaTilda_, radius_, regularisationPDE_, and regularise().

| void postProcessSens | ( | scalarField & | sens | ) |
Update part of fieldRegularisation to the sensitivitiy derivatives.
Definition at line 147 of file fieldRegularisation.C.
References betaArg_, mesh_, project_, regularise(), regularise_, and sharpenFunction_.

|
protected |
Definition at line 75 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), and postProcessSens().
|
protected |
Definition at line 76 of file fieldRegularisation.H.
Referenced by fieldRegularisation().
|
protected |
Cell zones related to topology optimisation.
Definition at line 81 of file fieldRegularisation.H.
Referenced by fieldRegularisation().
|
protected |
Perform regulaisation on alpha before inputing it on beta?
Definition at line 86 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), postProcessSens(), shouldRegularise(), and updateBeta().
|
protected |
Perform the projection (sharpening) step?
Definition at line 91 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), postProcessSens(), and updateBeta().
|
protected |
Smoothing radius.
May be isotropic or differ per spatial direction
Definition at line 98 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), and regularise().
|
protected |
Alpha field (design variables of topology optimisation).
Definition at line 103 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), and updateBeta().
|
protected |
The regularised alpha field, if regulatisation is performed.
Definition at line 108 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), regularise(), regularise(), and updateBeta().
|
protected |
Function used to sharpen the field after regularisation.
Definition at line 113 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), postProcessSens(), and updateBeta().
|
protected |
PDE used for the regularisation.
Definition at line 118 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), regularise(), and regularise().
|
protected |
Argument of the beta field.
Is either alpha_ (no regularisation) or alphaTilda_ (with regularisation)
Definition at line 126 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), postProcessSens(), and updateBeta().
|
protected |
Whether to apply a fixedValue BC or zeroGradient one to alphaTilda, when regularisation is performed.
Definition at line 132 of file fieldRegularisation.H.
Referenced by fieldRegularisation(), and growFromWalls().
|
protected |
Beta is the field used for all interpolations between fluid and solid in topology optimisation.
Definition at line 138 of file fieldRegularisation.H.
Referenced by beta(), fieldRegularisation(), and updateBeta().