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  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  (
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  {
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 
122 scalar 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 
150 void 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 
161 void objectiveFlowRatePartition::update_boundarydJdvn()
162 {
163  forAll(outletPatches_, pI)
164  {
165  const label patchI = outletPatches_[pI];
166  bdJdvnPtr_()[patchI] = flowRateDifference_[pI]/inletFlowRate_;
167  }
168 }
169 
170 
171 void 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];
180  objFunctionFilePtr_()
181  << setw(width_) << word("#" + patch.name() + "Tar") << " "
182  << setw(width_) << targetFlowRateFraction_[pI] << endl;
183  }
184 }
185 
186 
187 void objectiveFlowRatePartition::addHeaderColumns() const
188 {
189  for (const label patchI : outletPatches_)
190  {
191  const fvPatch& patch = mesh_.boundary()[patchI];
192  objFunctionFilePtr_()
193  << setw(width_) << word(patch.name() + "Ratio") << " ";
194  }
195 }
196 
197 
198 void 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 // ************************************************************************* //
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
bool readIfPresentCompat(const word &keyword, std::initializer_list< std::pair< const char *, int >> compat, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val using any compatibility names if needed. FatalIOError if it is found and there are excess tokens.
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:91
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
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:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
defineTypeNameAndDebug(objectivePartialVolume, 1)
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
Istream and Ostream manipulators taking arguments.
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Abstract base class for objective functions in incompressible flows.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127