57void objectivePowerDissipation::populateFieldNames()
66 adjointRAS().getAdjointTMVariablesBaseNames();
82 const word& adjointSolverName,
83 const word& primalSolverName
86 objectiveIncompressible(
mesh,
dict, adjointSolverName, primalSolverName),
87 zones_(mesh_.cellZones().indices(
dict.get<wordRes>(
"zones")))
94 extendedVariableName(
"Ua")
136 (
"gradDxdbMult" +
type()),
159 for (
const label zI : zones_)
165 J_ += 0.5*
gSum(integrandZone*VZone);
171 scalar porosityContr =
Zero;
172 for (
const label cellI : zoneI)
174 porosityContr +=
beta[cellI]*
magSqr(
U[cellI])*V[cellI];
196 const incompressibleAdjointSolver& adjSolver =
199 const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointRAS =
200 adjSolver.getAdjointVars().adjointTurbulence();
201 tmp<volScalarField> dnutdUMult = 0.5*
magSqr(S());
202 tmp<volVectorField> dnutdU = adjointRAS->nutJacobianU(dnutdUMult);
205 for (
const label zI : zones_)
208 for (
const label cellI : zoneI)
210 dJdvPtr_()[cellI] = dnutdU()[cellI];
218 for (
const label zI : zones_)
220 const cellZone& zoneI =
mesh_.cellZones()[zI];
221 for (
const label cellI : zoneI)
223 dJdvPtr_()[cellI] += integrand[cellI];
228 if (
mesh_.foundObject<topOVariablesBase>(
"topOVars"))
230 const topOVariablesBase& vars =
231 mesh_.lookupObject<topOVariablesBase>(
"topOVars");
233 const scalar betaMax = vars.getBetaMax();
234 for (
const label zI : zones_)
236 const cellZone& zoneI =
mesh_.cellZones()[zI];
237 for (
const label cellI : zoneI)
282 for (
const label zI : zones_)
285 for (
const label cellI : zoneI)
287 divDxDbMult[cellI] = integrand[cellI];
303 for (
const label zI : zones_)
306 for (
const label cellI : zoneI)
308 gradDxDbMult[cellI] = integrand[cellI];
311 gradDxDbMult.correctBoundaryConditions();
322 const topOVariablesBase& vars =
324 const scalar betaMax = vars.getBetaMax();
325 for (
const label zI : zones_)
328 for (
const label cellI : zoneI)
341 populateFieldNames();
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
compressible::turbulenceModel & turb
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
void correctBoundaryConditions()
Correct boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
void push_back(const T &val)
Append an element at the end of the list.
void setSize(label n)
Alias for resize().
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
void size(const label n)
Older name for setAddressableSize.
label find(const T &val) const
Find index of the first occurrence of the value.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Base class for selecting the betaMax value, i.e. the value multiplying the Brinkman penalisation term...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Mesh data needed to do the Finite Volume discretisation.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Time & time() const
Return the top-level database.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Base class for incompressibleAdjoint solvers.
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
const volVectorField & UInst() const
Return const reference to velocity.
const volVectorField & U() const
Return const reference to velocity.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Abstract base class for objective functions in incompressible flows.
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
const incompressibleVars & vars_
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
autoPtr< volVectorField > dJdvPtr_
const word objectiveName_
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
virtual scalar weight() const
Return the objective function weight.
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
const volScalarField & dJdb() const
Contribution to field sensitivities.
scalar J_
Objective function value and weight.
const dictionary & dict() const
Return objective dictionary.
const word adjointSolverName_
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Computes and minimizes the power dissipation within given cellZones. In the absence of significant vi...
virtual void update_dJdb()
Contribution to field sensitivities.
virtual void update_gradDxDbMultiplier()
Update grad(dx/db multiplier). Volume-based sensitivity term.
objectivePowerDissipation(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
virtual void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
virtual void addSource(fvScalarMatrix &matrix)
Add source terms to the adjoint turbulence model equations.
virtual scalar J()
Return the objective function value.
virtual void update_dJdTMvar1()
Update field to be added to the first adjoint turbulence model PDE.
virtual void update_divDxDbMultiplier()
Update div(dx/db multiplier). Volume-based sensitivity term.
virtual void update_dJdv()
Update values to be added to the adjoint outlet velocity.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
word extendedVariableName(const word &varName) const
Given a variable name, return a name that is possibly appended by the solverName, depending on useSol...
A class for managing temporary objects.
Base class for all design variables related to topology optimisation (topO). Provides the lookup func...
scalar getBetaMax() const
Get betaMax value.
virtual const volScalarField & beta() const =0
Get field used for physical interpolations.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
List< word > wordList
List of word.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
static constexpr const zero Zero
Global zero (0).
autoPtr< GeometricField< Type, fvPatchField, volMesh > > createZeroFieldPtr(const fvMesh &mesh, const word &name, const dimensionSet dims, bool printAllocation=false)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.