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 -------------------------------------------------------------------------------
12 License
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 
38 namespace Foam
39 {
40 
41 namespace objectives
42 {
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(objectivePtLosses, 1);
48 (
49  objectiveIncompressible,
50  objectivePtLosses,
51  dictionary
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());
101  const surfaceScalarField& phi = vars_.phiInst();
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  {
241  objFunctionFilePtr_()
242  << setw(width_) << mesh_.boundary()[patchI].name() << " ";
243  }
244 }
245 
246 
248 {
249  for (const scalar pt : patchPt_)
250  {
251  objFunctionFilePtr_()
252  << setw(width_) << pt << " ";
253  }
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 } // End namespace objectives
260 } // End namespace Foam
261 
262 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< vector >
Foam::objectives::objectivePtLosses::initialize
void initialize()
Return the objectiveReportPatches.
Definition: objectivePtLosses.C:82
Foam::incompressibleVars::phiInst
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
Definition: incompressibleVars.C:406
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::objectives::objectivePtLosses::addColumnValues
virtual void addColumnValues() const
Write information to additional columns.
Definition: objectivePtLosses.C:247
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< label >
objectivePtLosses.H
Foam::objectiveIncompressible::vars_
const incompressibleVars & vars_
Definition: objectiveIncompressible.H:62
Foam::objectives::addToRunTimeSelectionTable
addToRunTimeSelectionTable(objectiveIncompressible, objectiveForce, dictionary)
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::objectiveIncompressible::bdJdvnPtr_
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
Definition: objectiveIncompressible.H:82
Foam::objectives::objectivePtLosses::update_boundarydJdv
void update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
Definition: objectivePtLosses.C:180
Foam::objectiveIncompressible::bdJdvPtr_
autoPtr< boundaryVectorField > bdJdvPtr_
Definition: objectiveIncompressible.H:79
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::incompressibleVars::UInst
const volVectorField & UInst() const
Return const reference to velocity.
Definition: incompressibleVars.C:394
Foam::objectiveIncompressible::bdJdvtPtr_
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
Definition: objectiveIncompressible.H:85
Foam::incompressibleVars::p
const volScalarField & p() const
Return const reference to pressure.
Definition: incompressibleVars.C:305
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
coupledFvPatch.H
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::objectives::defineTypeNameAndDebug
defineTypeNameAndDebug(objectiveForce, 0)
Foam::Field< vector >
Foam::objective::J_
scalar J_
Objective function value and weight.
Definition: objective.H:77
createZeroField.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::incompressibleVars::U
const volVectorField & U() const
Return const reference to velocity.
Definition: incompressibleVars.C:331
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::objectiveIncompressible
Abstract base class for objective functions in incompressible flows.
Definition: objectiveIncompressible.H:54
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::objectives::objectivePtLosses::addHeaderColumns
virtual void addHeaderColumns() const
Write headers for additional columns.
Definition: objectivePtLosses.C:237
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
U
U
Definition: pEqn.H:72
Foam::objectives::objectivePtLosses::J
scalar J()
Return the objective function value.
Definition: objectivePtLosses.C:135
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::objectives::objectivePtLosses::update_boundarydJdp
void update_boundarydJdp()
Update values to be added to the adjoint inlet velocity.
Definition: objectivePtLosses.C:164
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::objective::mesh_
const fvMesh & mesh_
Definition: objective.H:67
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::polyBoundaryMesh::patchSet
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.
Definition: polyBoundaryMesh.C:864
Foam::objectiveIncompressible::bdJdpPtr_
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
Definition: objectiveIncompressible.H:88
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::objectives::objectivePtLosses::objectivePtLosses
objectivePtLosses(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
Definition: objectivePtLosses.C:58
Foam::objectives::objectivePtLosses::update_boundarydJdvn
void update_boundarydJdvn()
Update values to be added to the adjoint outlet pressure.
Definition: objectivePtLosses.C:200
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::objectives::objectivePtLosses::update_boundarydJdvt
void update_boundarydJdvt()
Update values to be added to the adjoint outlet tangential velocity.
Definition: objectivePtLosses.C:220
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::incompressibleVars::pInst
const volScalarField & pInst() const
Return const reference to pressure.
Definition: incompressibleVars.C:382