objectivePtLosses.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-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 FOSS GP
10 Copyright (C) 2019-2020 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "objectivePtLosses.H"
31#include "createZeroField.H"
32#include "coupledFvPatch.H"
33#include "IOmanip.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41namespace objectives
42{
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
48(
52);
53
54
55// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56
58(
59 const fvMesh& mesh,
60 const dictionary& dict,
61 const word& adjointSolverName,
62 const word& primalSolverName
63)
64:
65 objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
66 patches_(0),
67 patchPt_(0)
68{
69 // Find inlet/outlet patches
70 initialize();
71
72 // Allocate boundary field pointers
73 bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
74 bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
75 bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
76 bdJdvtPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
77}
78
79
80// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81
83{
84 // If patches are prescribed, use them
85
86 wordRes patchSelection;
87 if (dict().readIfPresent("patches", patchSelection))
88 {
89 patches_ = mesh_.boundaryMesh().patchSet(patchSelection).sortedToc();
90 }
91 // Otherwise, pick them up based on the mass flow.
92 // Note: a non-zero U initialisation should be used in order to pick up the
93 // outlet patches correctly
94 else
95 {
97 << "No patches provided to PtLosses. "
98 << "Choosing them according to the patch mass flows" << nl;
99
100 DynamicList<label> objectiveReportPatches(mesh_.boundary().size());
102 forAll(mesh_.boundary(), patchI)
103 {
104 const fvsPatchScalarField& phiPatch = phi.boundaryField()[patchI];
105 if (!isA<coupledFvPatch>(mesh_.boundary()[patchI]))
106 {
107 const scalar mass = gSum(phiPatch);
108 if (mag(mass) > SMALL)
109 {
110 objectiveReportPatches.append(patchI);
111 }
112 }
113 }
114 patches_.transfer(objectiveReportPatches);
115 }
116 patchPt_.setSize(patches_.size());
117
118 if (patches_.empty())
119 {
121 << "No valid patch name on which to minimize " << type() << endl
122 << exit(FatalError);
123 }
124 if (debug)
125 {
126 Info<< "Minimizing " << type() << " in patches:" << endl;
127 forAll(patches_, pI)
128 {
129 Info<< "\t " << mesh_.boundary()[patches_[pI]].name() << endl;
130 }
131 }
132}
133
134
136{
137 J_ = Zero;
138
139 // References
140 const volScalarField& p = vars_.pInst();
141 const volVectorField& U = vars_.UInst();
142
143 // Inlet/outlet patches
144 forAll(patches_, oI)
145 {
146 const label patchI = patches_[oI];
147 const vectorField& Sf = mesh_.boundary()[patchI].Sf();
148 scalar pt = -gSum
149 (
150 (U.boundaryField()[patchI] & Sf)
151 *(
152 p.boundaryField()[patchI]
153 + 0.5*magSqr(U.boundaryField()[patchI])
154 )
155 );
156 patchPt_[oI] = mag(pt);
157 J_ += pt;
158 }
159
160 return J_;
161}
162
163
165{
166 const volVectorField& U = vars_.U();
167
168 forAll(patches_, oI)
169 {
170 const label patchI = patches_[oI];
171
172 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
173 const vectorField& nf = tnf();
174
175 bdJdpPtr_()[patchI] = -(U.boundaryField()[patchI] & nf)*nf;
176 }
177}
178
179
181{
182 const volScalarField& p = vars_.p();
183 const volVectorField& U = vars_.U();
184
185 forAll(patches_, oI)
186 {
187 const label patchI = patches_[oI];
188
189 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
190 const vectorField& nf = tnf();
191 const fvPatchVectorField& Ub = U.boundaryField()[patchI];
192
193 bdJdvPtr_()[patchI] =
194 - (p.boundaryField()[patchI] + 0.5*magSqr(Ub))*nf
195 - (Ub & nf)*Ub;
196 }
197}
198
199
201{
202 const volScalarField& p = vars_.p();
203 const volVectorField& U = vars_.U();
204
205 forAll(patches_, oI)
206 {
207 const label patchI = patches_[oI];
208
209 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
210 const vectorField& nf = tnf();
211
212 bdJdvnPtr_()[patchI] =
213 - p.boundaryField()[patchI]
214 - 0.5*magSqr(U.boundaryField()[patchI])
215 - sqr(U.boundaryField()[patchI] & nf);
216 }
217}
218
219
221{
222 const volVectorField& U = vars_.U();
223
224 forAll(patches_, oI)
225 {
226 const label patchI = patches_[oI];
227
228 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
229 const vectorField& nf = tnf();
230 scalarField Un(U.boundaryField()[patchI] & nf);
231
232 bdJdvtPtr_()[patchI] = -Un*(U.boundaryField()[patchI] - Un*nf);
233 }
234}
235
236
238{
239 for (const label patchI : patches_)
240 {
242 << setw(width_) << mesh_.boundary()[patchI].name() << " ";
243 }
244}
245
246
248{
249 for (const scalar pt : patchPt_)
250 {
252 << setw(width_) << pt << " ";
253 }
254}
255
256
257// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258
259} // End namespace objectives
260} // End namespace Foam
261
262// ************************************************************************* //
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.
surfaceScalarField & phi
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
const Boundary & boundaryField() const
Return const-reference to the boundary field.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const volVectorField & UInst() const
Return const reference to velocity.
const volVectorField & U() const
Return const reference to velocity.
const volScalarField & p() const
Return const reference to pressure.
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
const volScalarField & pInst() const
Return const reference to pressure.
Abstract base class for objective functions in incompressible flows.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
const incompressibleVars & vars_
autoPtr< boundaryVectorField > bdJdvPtr_
const fvMesh & mesh_
Definition: objective.H:67
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:141
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:152
scalar J_
Objective function value and weight.
Definition: objective.H:77
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:94
void initialize()
Return the objectiveReportPatches.
virtual void addHeaderColumns() const
Write headers for additional columns.
void update_boundarydJdvn()
Update values to be added to the adjoint outlet pressure.
scalar J()
Return the objective function value.
void update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
virtual void addColumnValues() const
Write information to additional columns.
void update_boundarydJdp()
Update values to be added to the adjoint inlet velocity.
void update_boundarydJdvt()
Update values to be added to the adjoint outlet tangential velocity.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A class for managing temporary objects.
Definition: tmp.H:65
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333