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 -------------------------------------------------------------------------------
11 License
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 
36 namespace Foam
37 {
38 
39 namespace objectives
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(objectiveFlowRatePartition, 0);
46 (
47  objectiveIncompressible,
48  objectiveFlowRatePartition,
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  targetFlowRatePercentage_(),
79  currentFlowRatePercentage_(outletPatches_.size(), Zero),
80  inletFlowRate_(0),
81  flowRateDifference_(outletPatches_.size(), Zero)
82 {
83  // Read target percentages if present, otherwise treat them as equally
84  // distributed
85  if (!dict.readIfPresent("targetPercentages", targetFlowRatePercentage_))
86  {
87  const label nOutPatches = outletPatches_.size();
88  targetFlowRatePercentage_.setSize(nOutPatches, 1.0/scalar(nOutPatches));
89  }
90  // Sanity checks
91  if (targetFlowRatePercentage_.size() != outletPatches_.size())
92  {
94  << "Inconsistent sizes for targetPercentages and outletPatches"
95  << exit(FatalError);
96 
97  }
98 
99  // Allocate boundary field pointers
100  bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
101  bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
102 
103  // Determine word width for output
104  for (const label patchI : outletPatches_)
105  {
106  const fvPatch& patch = mesh_.boundary()[patchI];
107  unsigned int wordSize = word(patch.name() + "Ratio").size();
108  width_ = max(width_, wordSize);
109  }
110 }
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 {
117  J_ = 0;
118  const volVectorField& U = vars_.UInst();
119 
120  // Inlet patches
121  inletFlowRate_ = 0;
122  for (const label patchI : inletPatches_)
123  {
124  const vectorField& Sf = mesh_.boundary()[patchI].Sf();
125  // Negative value
126  inletFlowRate_ += gSum(U.boundaryField()[patchI] & Sf);
127  }
128 
129  // Outlet patches
130  forAll(outletPatches_, pI)
131  {
132  const label patchI = outletPatches_[pI];
133  const vectorField& Sf = mesh_.boundary()[patchI].Sf();
134  const scalar outletFlowRate = gSum(U.boundaryField()[patchI] & Sf);
135  currentFlowRatePercentage_[pI] = -outletFlowRate/inletFlowRate_;
136  flowRateDifference_[pI] =
137  targetFlowRatePercentage_[pI] - currentFlowRatePercentage_[pI];
138  J_ += 0.5*flowRateDifference_[pI]*flowRateDifference_[pI];
139  }
140 
141  return J_;
142 }
143 
144 
146 {
147  forAll(outletPatches_, pI)
148  {
149  const label patchI = outletPatches_[pI];
150  tmp<vectorField> nf = mesh_.boundary()[patchI].nf();
151  bdJdvPtr_()[patchI] = nf*flowRateDifference_[pI]/inletFlowRate_;
152  }
153 }
154 
155 
157 {
158  forAll(outletPatches_, pI)
159  {
160  const label patchI = outletPatches_[pI];
161  bdJdvnPtr_()[patchI] = flowRateDifference_[pI]/inletFlowRate_;
162  }
163 }
164 
165 
167 {
169  << setw(width_) << "#inletFlowRate" << " "
170  << setw(width_) << -inletFlowRate_ << endl;
171  forAll(outletPatches_, pI)
172  {
173  const label patchI = outletPatches_[pI];
174  const fvPatch& patch = mesh_.boundary()[patchI];
176  << setw(width_) << word("#" + patch.name() + "Tar") << " "
177  << setw(width_) << targetFlowRatePercentage_[pI] << endl;
178  }
179 }
180 
181 
183 {
184  for (const label patchI : outletPatches_)
185  {
186  const fvPatch& patch = mesh_.boundary()[patchI];
188  << setw(width_) << word(patch.name() + "Ratio") << " ";
189  }
190 }
191 
192 
194 {
195  for (const scalar flowRate : currentFlowRatePercentage_)
196  {
198  << setw(width_) << flowRate << " ";
199  }
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 } // End namespace objectives
206 } // End namespace Foam
207 
208 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
void update_boundarydJdvn()
Update values to be added to the adjoint outlet pressure.
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:87
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:173
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
objectiveFlowRatePartition(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
defineTypeNameAndDebug(objectiveFlowRate, 0)
virtual void addHeaderColumns() const
Write headers for additional columns.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
const fvMesh & mesh_
Definition: objective.H:63
scalar J()
Return the objective function value.
scalar J_
Objective function value and weight.
Definition: objective.H:75
autoPtr< boundaryVectorField > bdJdvPtr_
Istream and Ostream manipulators taking arguments.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
const incompressibleVars & vars_
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:190
addToRunTimeSelectionTable(objectiveIncompressible, objectiveFlowRate, dictionary)
virtual void addColumnValues() const
Write information to additional columns.
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const std::string patch
OpenFOAM patch number as a std::string.
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
Field< vector > vectorField
Specialisation of Field<T> for vector.
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const volVectorField & UInst() const
Return const reference to velocity.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
Abstract base class for objective functions in incompressible flows.
Namespace for OpenFOAM.
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:705
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157