A regulatisation PDE based on a Helmholtz filter. More...
#include <Helmholtz.H>


Public Member Functions | |
| TypeName ("Helmholtz") | |
| Runtime type information. | |
| Helmholtz (const fvMesh &mesh, const dictionary &dict, const topOZones &zones) | |
| Construct from components. | |
| virtual | ~Helmholtz ()=default |
| Destructor. | |
| virtual void | regularise (const volScalarField &aTilda, const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius, const scalar minSetValue=Zero, const bool fixATildaValues=true) |
| Regularise field (or sensitivities) using a Laplace PDE. | |
| Public Member Functions inherited from regularisationPDE | |
| TypeName ("regularisationPDE") | |
| Runtime type information. | |
| declareRunTimeSelectionTable (autoPtr, regularisationPDE, dictionary,(const fvMesh &mesh, const dictionary &dict, const topOZones &zones),(mesh, dict, zones)) | |
| regularisationPDE (const fvMesh &mesh, const dictionary &dict, const topOZones &zones) | |
| Construct from components. | |
| virtual | ~regularisationPDE ()=default |
| Destructor. | |
Protected Member Functions | |
| void | solveEqn (const volScalarField &aTilda, const scalarField &source, scalarField &result, const bool isTopoField, const regularisationRadius &radius, const scalar minSetValue=Zero, const bool fixATildaValues=true) |
| Solve the regulatisation equation. | |
| Protected Member Functions inherited from regularisationPDE | |
| void | setValues (fvScalarMatrix &bTildaEqn, const bool isTopoField, const scalar minSetValue=Zero) |
| Set fixed bTilda values. | |
| void | setValues (const fvMesh &mesh, DynamicList< label > &cells, DynamicList< scalar > &values, bool isTopoField, const scalar minSetValue=Zero) |
| Gather the cells with constant values from all constrained zones. | |
| scalar | computeRadius () |
| Compute smoothing radius, if not directly given. | |
Protected Attributes | |
| bool | solveOnActiveCells_ |
| Solve the regularisationPDE only on a subset mesh made of the active cell zones. | |
| scalar | wallValue_ |
| Fixed value at wall boundaries. Defaults to 1. | |
| Protected Attributes inherited from regularisationPDE | |
| const fvMesh & | mesh_ |
| const dictionary | dict_ |
| const topOZones & | zones_ |
| Cell zones related to topology optimisation. | |
| bool | growFromWalls_ |
| Whether to apply a fixedValue BC or zeroGradient one to alphaTilda, when regularisation is performed. | |
Additional Inherited Members | |
| Static Public Member Functions inherited from regularisationPDE | |
| static autoPtr< regularisationPDE > | New (const fvMesh &mesh, const dictionary &dict, const topOZones &zones) |
| Construct and return the selected regularisationPDE. | |
A regulatisation PDE based on a Helmholtz filter.
Definition at line 49 of file Helmholtz.H.
| Helmholtz | ( | const fvMesh & | mesh, |
| const dictionary & | dict, | ||
| const topOZones & | zones ) |
Construct from components.
Definition at line 150 of file Helmholtz.C.
References dict, mesh, solveOnActiveCells_, and wallValue_.
|
virtualdefault |
Destructor.
References Foam::Zero.
|
protected |
Solve the regulatisation equation.
The mesh used is the one of aTilda
Definition at line 43 of file Helmholtz.C.
References regularisationRadius::addRegularisationTerm(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), regularisationPDE::dict_, Foam::dimless, forAll, regularisationPDE::growFromWalls_, Foam::Info, SolverPerformance< Type >::initialResidual(), GeometricField< Type, PatchField, GeoMesh >::internalField(), Foam::isA(), iters(), Foam::mag(), DimensionedField< Type, GeoMesh >::mesh(), mesh, IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, GeometricField< Type, PatchField, GeoMesh >::primitiveField(), regularisationPDE::setValues(), fvMatrix< Type >::solve(), Foam::fvm::Sp(), wallValue_, and Foam::Zero.
Referenced by regularise().


|
virtual |
Regularise field (or sensitivities) using a Laplace PDE.
Implements regularisationPDE.
Definition at line 165 of file Helmholtz.C.
References cells, regularisationPDE::mesh_, Field< Type >::rmap(), regularisationPDE::setValues(), and solveEqn().

|
protected |
Solve the regularisationPDE only on a subset mesh made of the active cell zones.
Definition at line 77 of file Helmholtz.H.
Referenced by Helmholtz().
|
protected |
Fixed value at wall boundaries. Defaults to 1.
Definition at line 82 of file Helmholtz.H.
Referenced by Helmholtz(), and solveEqn().