59 const word& adjointSolverName,
60 const word& primalSolverName
63 objectiveIncompressible(
mesh,
dict, adjointSolverName, primalSolverName),
78 targetFlowRateFraction_(),
79 currentFlowRateFraction_(outletPatches_.
size(),
Zero),
81 flowRateDifference_(outletPatches_.
size(),
Zero)
87 !
dict.readIfPresentCompat
89 "targetFractions", {{
"targetPercentages", 2306}},
90 targetFlowRateFraction_
94 const label nOutPatches = outletPatches_.size();
95 targetFlowRateFraction_.setSize(nOutPatches, 1.0/scalar(nOutPatches));
98 if (targetFlowRateFraction_.size() != outletPatches_.size())
101 <<
"Inconsistent sizes for targetFractions and outletPatches"
111 for (
const label patchI : outletPatches_)
113 const fvPatch& patch = mesh_.boundary()[patchI];
114 unsigned int wordSize =
word(patch.name() +
"Ratio").size();
122scalar objectiveFlowRatePartition::J()
129 for (
const label patchI : inletPatches_)
132 inletFlowRate_ +=
gSum(
phi.boundaryField()[patchI]);
136 forAll(outletPatches_, pI)
138 const label patchI = outletPatches_[pI];
139 const scalar outletFlowRate =
gSum(
phi.boundaryField()[patchI]);
140 currentFlowRateFraction_[pI] = -outletFlowRate/inletFlowRate_;
141 flowRateDifference_[pI] =
142 targetFlowRateFraction_[pI] - currentFlowRateFraction_[pI];
143 J_ += 0.5*flowRateDifference_[pI]*flowRateDifference_[pI];
150void objectiveFlowRatePartition::update_boundarydJdv()
152 forAll(outletPatches_, pI)
154 const label patchI = outletPatches_[pI];
156 bdJdvPtr_()[patchI] = nf*flowRateDifference_[pI]/inletFlowRate_;
161void objectiveFlowRatePartition::update_boundarydJdvn()
163 forAll(outletPatches_, pI)
165 const label patchI = outletPatches_[pI];
166 bdJdvnPtr_()[patchI] = flowRateDifference_[pI]/inletFlowRate_;
171void objectiveFlowRatePartition::addHeaderInfo()
const
173 objFunctionFilePtr_()
174 <<
setw(width_) <<
"#inletFlowRate" <<
" "
175 <<
setw(width_) << -inletFlowRate_ <<
endl;
176 forAll(outletPatches_, pI)
178 const label patchI = outletPatches_[pI];
179 const fvPatch& patch = mesh_.boundary()[patchI];
187void objectiveFlowRatePartition::addHeaderColumns()
const
189 for (
const label patchI : outletPatches_)
198void objectiveFlowRatePartition::addColumnValues()
const
200 for (
const scalar flowRate : currentFlowRateFraction_)
202 objFunctionFilePtr_()
203 <<
setw(width_) << flowRate <<
" ";
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.
for(const label curEdgei :curPointEdges)
label size() const noexcept
The number of elements in list.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
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.
wordList sortedToc() const
Return the sorted table of contents.
Mesh data needed to do the Finite Volume discretisation.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Abstract base class for objective functions in incompressible flows.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
const incompressibleVars & vars_
autoPtr< boundaryVectorField > bdJdvPtr_
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
unsigned int width_
Default width of entries when writing in the objective files.
scalar J_
Objective function value and weight.
const dictionary & dict() const
Return objective dictionary.
Split inlet flow rate to given percentages at the prescribed outlet patches.
objectiveFlowRatePartition(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
virtual void addHeaderColumns() const
Write headers for additional columns.
virtual void update_boundarydJdvn()
Update values to be added to the adjoint outlet pressure.
virtual scalar J()
Return the objective function value.
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
virtual void update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
virtual void addColumnValues() const
Write information to additional columns.
A class for managing temporary objects.
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.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Omanip< int > setw(const int i)
Ostream & endl(Ostream &os)
Add newline and flush stream.
static constexpr const zero Zero
Global zero (0).
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
#define forAll(list, i)
Loop across all elements in list.