objectiveFlowRate.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 
29 #include "objectiveFlowRate.H"
30 #include "createZeroField.H"
31 #include "IOmanip.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace objectives
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(objectiveFlowRate, 0);
45 (
46  objectiveIncompressible,
47  objectiveFlowRate,
49 );
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const fvMesh& mesh,
57  const dictionary& dict,
58  const word& adjointSolverName,
59  const word& primalSolverName
60 )
61 :
62  objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
63  patches_
64  (
65  mesh_.boundaryMesh().patchSet(dict.get<wordRes>("patches")).sortedToc()
66  ),
67  flowRates_(patches_.size(), Zero)
68 {
69  // Allocate boundary field pointers
70  bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
71  bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
72 }
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
77 scalar objectiveFlowRate::J()
78 {
79  J_ = 0;
81 
82  forAll(patches_, pI)
83  {
84  const label patchI = patches_[pI];
85  flowRates_[pI] = gSum(phi.boundaryField()[patchI]);
86  J_ += flowRates_[pI];
87  }
88 
89  return J_;
90 }
91 
92 
94 {
95  for (const label patchI : patches_)
96  {
97  bdJdvPtr_()[patchI] = mesh_.boundary()[patchI].nf();
98  }
99 }
100 
101 
103 {
104  for (const label patchI : patches_)
105  {
106  bdJdvnPtr_()[patchI] = 1;
107  }
108 }
109 
110 
112 {
113  for (const label patchI : patches_)
114  {
115  const fvPatch& patch = mesh_.boundary()[patchI];
117  << setw(width_) << word(patch.name() + "FlowRate") << " ";
118  }
119 }
120 
121 
123 {
124  for (const scalar flowRate : flowRates_)
125  {
127  << setw(width_) << flowRate << " ";
128  }
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 } // End namespace objectives
135 } // End namespace Foam
136 
137 // ************************************************************************* //
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
objectiveFlowRate(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:209
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
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:421
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
virtual scalar J()
Return the objective function value.
virtual void update_boundarydJdvn()
Update values to be added to the adjoint outlet pressure.
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
const fvMesh & mesh_
Definition: objective.H:63
scalar J_
Objective function value and weight.
Definition: objective.H:76
autoPtr< boundaryVectorField > bdJdvPtr_
Istream and Ostream manipulators taking arguments.
virtual void addColumnValues() const
Write information to additional columns.
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
const incompressibleVars & vars_
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:226
virtual void update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
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
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
Abstract base class for objective functions in incompressible flows.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127