Loading...
Searching...
No Matches
objectiveFlowRatePartition.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2022 PCOpt/NTUA
9 Copyright (C) 2013-2022 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "createZeroField.H"
31#include "IOmanip.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
39namespace objectives
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
46(
50);
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const fvMesh& mesh,
58 const dictionary& dict,
59 const word& adjointSolverName,
60 const word& primalSolverName
61)
62:
63 objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
64 inletPatches_
65 (
66 mesh_.boundaryMesh().patchSet
67 (
68 dict.get<wordRes>("inletPatches")
69 ).sortedToc()
70 ),
71 outletPatches_
72 (
73 mesh_.boundaryMesh().patchSet
74 (
75 dict.get<wordRes>("outletPatches")
76 ).sortedToc()
77 ),
78 targetFlowRateFraction_(),
79 currentFlowRateFraction_(outletPatches_.size(), Zero),
80 inletFlowRate_(0),
81 flowRateDifference_(outletPatches_.size(), Zero)
82{
83 // Read target fractions if present, otherwise treat them as equally
84 // distributed
85 if
86 (
87 !dict.readIfPresentCompat
88 (
89 "targetFractions", {{"targetPercentages", 2306}},
90 targetFlowRateFraction_
91 )
92 )
93 {
94 const label nOutPatches = outletPatches_.size();
95 targetFlowRateFraction_.setSize(nOutPatches, 1.0/scalar(nOutPatches));
96 }
97 // Sanity checks
98 if (targetFlowRateFraction_.size() != outletPatches_.size())
99 {
100 FatalErrorInFunction
101 << "Inconsistent sizes for targetFractions and outletPatches"
102 << exit(FatalError);
103
104 }
105
106 // Allocate boundary field pointers
107 bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
108 bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
109
110 // Determine word width for output
111 for (const label patchI : outletPatches_)
112 {
113 const fvPatch& patch = mesh_.boundary()[patchI];
114 unsigned int wordSize = word(patch.name() + "Ratio").size();
115 width_ = max(width_, wordSize);
116 }
117}
118
119
120// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121
122scalar objectiveFlowRatePartition::J()
123{
124 J_ = 0;
125 const surfaceScalarField& phi = vars_.phiInst();
126
127 // Inlet patches
128 inletFlowRate_ = 0;
129 for (const label patchI : inletPatches_)
130 {
131 // Negative value
132 inletFlowRate_ += gSum(phi.boundaryField()[patchI]);
133 }
134
135 // Outlet patches
136 forAll(outletPatches_, pI)
137 {
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];
144 }
145
146 return J_;
147}
148
149
150void objectiveFlowRatePartition::update_boundarydJdv()
151{
152 forAll(outletPatches_, pI)
153 {
154 const label patchI = outletPatches_[pI];
155 tmp<vectorField> nf = mesh_.boundary()[patchI].nf();
156 bdJdvPtr_()[patchI] = nf*flowRateDifference_[pI]/inletFlowRate_;
157 }
158}
159
160
161void objectiveFlowRatePartition::update_boundarydJdvn()
162{
163 forAll(outletPatches_, pI)
165 const label patchI = outletPatches_[pI];
166 bdJdvnPtr_()[patchI] = flowRateDifference_[pI]/inletFlowRate_;
167 }
168}
169
170
171void objectiveFlowRatePartition::addHeaderInfo() const
172{
173 objFunctionFilePtr_()
174 << setw(width_) << "#inletFlowRate" << " "
175 << setw(width_) << -inletFlowRate_ << endl;
176 forAll(outletPatches_, pI)
177 {
178 const label patchI = outletPatches_[pI];
179 const fvPatch& patch = mesh_.boundary()[patchI];
181 << setw(width_) << word("#" + patch.name() + "Tar") << " "
182 << setw(width_) << targetFlowRateFraction_[pI] << endl;
183 }
184}
185
186
187void objectiveFlowRatePartition::addHeaderColumns() const
188{
189 for (const label patchI : outletPatches_)
190 {
191 const fvPatch& patch = mesh_.boundary()[patchI];
193 << setw(width_) << word(patch.name() + "Ratio") << " ";
194 }
195}
196
197
198void objectiveFlowRatePartition::addColumnValues() const
199{
200 for (const scalar flowRate : currentFlowRateFraction_)
201 {
202 objFunctionFilePtr_()
203 << setw(width_) << flowRate << " ";
204 }
205}
206
207
208// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209
210} // End namespace objectives
211} // End namespace Foam
212
213// ************************************************************************* //
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.
Definition DLListBase.H:194
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,...
Definition dictionary.H:133
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.
Definition dictionary.C:601
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Abstract base class for objective functions in incompressible flows.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
const incompressibleVars & vars_
autoPtr< boundaryVectorField > bdJdvPtr_
const fvMesh & mesh_
Definition objective.H:63
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition objective.H:209
unsigned int width_
Default width of entries when writing in the objective files.
Definition objective.H:226
scalar J_
Objective function value and weight.
Definition objective.H:76
const dictionary & dict() const
Return objective dictionary.
Definition objective.C:91
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.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
Type gSum(const FieldField< Field, Type > &f)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299