Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
29#include "objectiveFlowRate.H"
30#include "createZeroField.H"
31#include "IOmanip.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
38namespace objectives
39{
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
45(
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
72}
73
74
75// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76
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// ************************************************************************* //
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.
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
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
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
Minimize/maximize flow rate through a given set of patches.
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 update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
virtual void addColumnValues() const
Write information to additional columns.
objectiveFlowRate(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
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.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Omanip< int > setw(const int i)
Definition IOmanip.H:199
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