Loading...
Searching...
No Matches
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
41namespace objectives
42{
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
48(
52);
53
54
55// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56
57void objectivePtLosses::initialize()
58{
59 // If patches are prescribed, use them
60
61 wordRes patchSelection;
62 if (dict().readIfPresent("patches", patchSelection))
63 {
64 patches_ = mesh_.boundaryMesh().patchSet(patchSelection).sortedToc();
65 }
66 // Otherwise, pick them up based on the mass flow.
67 // Note: a non-zero U initialisation should be used in order to pick up the
68 // outlet patches correctly
69 else
70 {
72 << "No patches provided to PtLosses. "
73 << "Choosing them according to the patch mass flows" << nl;
74
75 DynamicList<label> objectiveReportPatches(mesh_.boundary().size());
77 forAll(mesh_.boundary(), patchI)
78 {
79 const fvsPatchScalarField& phiPatch = phi.boundaryField()[patchI];
80 if (!isA<coupledFvPatch>(mesh_.boundary()[patchI]))
81 {
82 const scalar mass = gSum(phiPatch);
83 if (mag(mass) > SMALL)
84 {
85 objectiveReportPatches.push_back(patchI);
86 }
87 }
88 }
89 patches_.transfer(objectiveReportPatches);
90 }
91 patchPt_.setSize(patches_.size());
92
93 if (patches_.empty())
94 {
96 << "No valid patch name on which to minimize " << type() << endl
97 << exit(FatalError);
98 }
99 if (debug)
100 {
101 Info<< "Minimizing " << type() << " in patches:" << endl;
102 forAll(patches_, pI)
103 {
104 Info<< "\t " << mesh_.boundary()[patches_[pI]].name() << endl;
106 }
107}
108
109
110// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111
113(
114 const fvMesh& mesh,
115 const dictionary& dict,
116 const word& adjointSolverName,
117 const word& primalSolverName
118)
119:
120 objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
121 patches_(0),
122 patchPt_(0)
123{
124 // Find inlet/outlet patches
125 initialize();
126
127 // Allocate boundary field pointers
132}
133
134
135// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136
138{
139 J_ = Zero;
140
141 // References
142 const volScalarField& p = vars_.pInst();
143 const volVectorField& U = vars_.UInst();
144
145 // Inlet/outlet patches
146 forAll(patches_, oI)
147 {
148 const label patchI = patches_[oI];
149 const vectorField& Sf = mesh_.boundary()[patchI].Sf();
150 scalar pt = -gSum
151 (
152 (U.boundaryField()[patchI] & Sf)
153 *(
154 p.boundaryField()[patchI]
155 + 0.5*magSqr(U.boundaryField()[patchI])
156 )
157 );
158 patchPt_[oI] = mag(pt);
159 J_ += pt;
160 }
161
162 return J_;
163}
164
165
167{
168 const volVectorField& U = vars_.U();
169
170 forAll(patches_, oI)
171 {
172 const label patchI = patches_[oI];
173
174 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
175 const vectorField& nf = tnf();
176
177 bdJdpPtr_()[patchI] = -(U.boundaryField()[patchI] & nf)*nf;
178 }
179}
180
181
183{
184 const volScalarField& p = vars_.p();
185 const volVectorField& U = vars_.U();
186
187 forAll(patches_, oI)
188 {
189 const label patchI = patches_[oI];
190
191 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
192 const vectorField& nf = tnf();
193 const fvPatchVectorField& Ub = U.boundaryField()[patchI];
194
195 bdJdvPtr_()[patchI] =
196 - (p.boundaryField()[patchI] + 0.5*magSqr(Ub))*nf
197 - (Ub & nf)*Ub;
198 }
199}
200
201
203{
204 const volScalarField& p = vars_.p();
205 const volVectorField& U = vars_.U();
206
207 forAll(patches_, oI)
208 {
209 const label patchI = patches_[oI];
210
211 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
212 const vectorField& nf = tnf();
213
214 bdJdvnPtr_()[patchI] =
215 - p.boundaryField()[patchI]
216 - 0.5*magSqr(U.boundaryField()[patchI])
217 - sqr(U.boundaryField()[patchI] & nf);
218 }
219}
220
221
223{
224 const volVectorField& U = vars_.U();
225
226 forAll(patches_, oI)
227 {
228 const label patchI = patches_[oI];
229
230 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
231 const vectorField& nf = tnf();
232 scalarField Un(U.boundaryField()[patchI] & nf);
233
234 bdJdvtPtr_()[patchI] = -Un*(U.boundaryField()[patchI] - Un*nf);
235 }
236}
237
238
240{
241 for (const label patchI : patches_)
244 << setw(width_) << mesh_.boundary()[patchI].name() << " ";
245 }
246}
247
248
250{
251 for (const scalar pt : patchPt_)
252 {
254 << setw(width_) << pt << " ";
255 }
256}
257
258
259// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260
261} // End namespace objectives
262} // End namespace Foam
263
264// ************************************************************************* //
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.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void push_back(const T &val)
Copy append an element to the end of this list.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries 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:133
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...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
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: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
objectivePtLosses(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
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.
virtual void update_boundarydJdp()
Update values to be added to the adjoint inlet velocity.
virtual void update_boundarydJdvt()
Update values to be added to the adjoint outlet tangential velocity.
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A class for managing temporary objects.
Definition tmp.H:75
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
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< vector > fvPatchVectorField
fvsPatchField< scalar > fvsPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299