57#ifndef ShapeSensitivitiesBase_H
58#define ShapeSensitivitiesBase_H
80class ShapeSensitivitiesBase
82 public adjointSensitivity
184 ShapeSensitivitiesBase(
const ShapeSensitivitiesBase&) =
delete;
187 void operator=(
const ShapeSensitivitiesBase&) =
delete;
199 ShapeSensitivitiesBase
Useful typenames for fields defined only at the boundaries.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of boundary fields.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
tmp< volVectorField > getWallFaceSensVec()
Get wall face sensitivity vectors field.
autoPtr< pointBoundaryVectorField > wallPointSensNormalVecPtr_
Normal sens as vectors.
void clearSensitivities()
Zero sensitivity fields and their constituents.
bool hasMultiplier(bool(objective::*hasFunction)() const)
Check if any of the available objective has a certain multiplier, provided through a function object.
bool writeAllSurfaceFiles_
Whether to write all surface sensitivity fields.
tmp< volVectorField > getWallFaceSensNormalVec()
Get wall face normal sens as vectors field.
autoPtr< boundaryVectorField > wallFaceSensVecPtr_
Wall face sens w.r.t. (x,y.z).
labelHashSet sensitivityPatchIDs_
Patches on which to compute shape sensitivities.
virtual const boundaryScalarField & getWallFaceSensNormalBoundary() const
Get wall face sensitivity projected to normal field.
virtual const boundaryVectorField & getWallFaceSensNormalVecBoundary() const
Get wall face normal sens as vectors field.
void writeFaceBasedSens() const
Write face-based sensitivities, if present.
void writePointBasedSens() const
Write point-based sensitivities, if present.
virtual const boundaryVectorField & getWallFaceSensVecBoundary() const
Get wall face sensitivity vectors field.
TypeName("ShapeSensitivitiesBase")
Runtime type information.
autoPtr< pointBoundaryVectorField > wallPointSensVecPtr_
Wall point sens w.r.t. (x,y.z).
tmp< pointVectorField > getWallPointSensVec()
Get wall point sensitivity vectors field.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
tmp< pointScalarField > getWallPointSensNormal()
Get wall point sensitivity projected to normal field.
void allocateEikonalSolver()
Allocate the adjoint eikonal solver.
void clearSurfaceFields()
Clear surface/point fields.
virtual ~ShapeSensitivitiesBase()=default
Destructor.
tmp< GeometricField< Type, fvPatchField, volMesh > > constructVolSensitivtyField(const autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > &sensFieldPtr, const word &name) const
Constructs volField based on boundaryField and writes it.
void allocateMultipliers()
Allocate multiplier fields.
tmp< pointVectorField > getWallPointSensNormalVec()
Get wall point sens as vectors field.
autoPtr< pointBoundaryScalarField > wallPointSensNormalPtr_
Wall point sens projected to normal.
tmp< volScalarField > getWallFaceSensNormal()
Get wall face sensitivity projected to normal field.
void constructAndWriteSensitivityField(const autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > &sensFieldPtr, const word &name) const
Constructs volField based on boundaryField and writes it.
const labelHashSet & sensitivityPatchIDs() const
Get patch IDs on which sensitivities are computed.
autoPtr< boundaryScalarField > wallFaceSensNormalPtr_
Wall face sens projected to normal.
void setSensitivityPatchIDs(const labelHashSet &sensPatchIDs)
Overwrite sensitivityPatchIDs.
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
void constructAndWriteSensitivtyPointField(const autoPtr< List< Field< Type > > > &sensFieldPtr, const word &name) const
Constructs pointField based on boundaryField and writes it.
autoPtr< boundaryVectorField > wallFaceSensNormalVecPtr_
Normal sens as vectors.
virtual const labelHashSet & geometryVariationIntegrationPatches() const
Return set of patches on which to compute direct sensitivities.
void clearMultipliers()
Clear multipliers.
Base class for adjoint solvers.
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.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
const fvMesh & mesh() const
Return reference to mesh.
const dictionary & dict() const
Return the construction dictionary.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
static const word null
An empty word.
volScalarField::Boundary boundaryScalarField
volFields
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
volVectorField::Boundary boundaryVectorField
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.