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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 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  if (dict().found("patches"))
86  {
88  (
90  (
91  wordReList(dict().get<wordRes>("patches"))
92  )
93  );
94  patches_ = patches.toc();
95  }
96  // Otherwise, pick them up based on the mass flow.
97  // Note: a non-zero U initialisation should be used in order to pick up the
98  // outlet patches correctly
99  else
100  {
102  << "No patches provided to PtLosses. Chossing them according to "
103  << "the patch mass flows"
104  << endl;
105  DynamicList<label> objectiveReportPatches(mesh_.boundary().size());
106  const surfaceScalarField& phi = vars_.phiInst();
107  forAll(mesh_.boundary(), patchI)
108  {
109  const fvsPatchScalarField& phiPatch = phi.boundaryField()[patchI];
110  if (!isA<coupledFvPatch>(mesh_.boundary()[patchI]))
111  {
112  const scalar mass = gSum(phiPatch);
113  if (mag(mass) > SMALL)
114  {
115  objectiveReportPatches.append(patchI);
116  }
117  }
118  }
119  patches_.transfer(objectiveReportPatches);
120  }
121  patchPt_.setSize(patches_.size());
122 
123  if (patches_.empty())
124  {
126  << "No valid patch name on which to minimize " << type() << endl
127  << exit(FatalError);
128  }
129  if (debug)
130  {
131  Info<< "Minimizing " << type() << " in patches:" << endl;
132  forAll(patches_, pI)
133  {
134  Info<< "\t " << mesh_.boundary()[patches_[pI]].name() << endl;
135  }
136  }
137 }
138 
139 
141 {
142  J_ = Zero;
143 
144  // References
145  const volScalarField& p = vars_.pInst();
146  const volVectorField& U = vars_.UInst();
147 
148  // Inlet/outlet patches
149  forAll(patches_, oI)
150  {
151  const label patchI = patches_[oI];
152  const vectorField& Sf = mesh_.boundary()[patchI].Sf();
153  scalar pt = -gSum
154  (
155  (U.boundaryField()[patchI] & Sf)
156  *(
157  p.boundaryField()[patchI]
158  + 0.5*magSqr(U.boundaryField()[patchI])
159  )
160  );
161  patchPt_[oI] = mag(pt);
162  J_ += pt;
163  }
164 
165  return J_;
166 }
167 
168 
170 {
171  const volVectorField& U = vars_.U();
172 
173  forAll(patches_, oI)
174  {
175  const label patchI = patches_[oI];
176 
177  tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
178  const vectorField& nf = tnf();
179 
180  bdJdpPtr_()[patchI] = -(U.boundaryField()[patchI] & nf)*nf;
181  }
182 }
183 
184 
186 {
187  const volScalarField& p = vars_.p();
188  const volVectorField& U = vars_.U();
189 
190  forAll(patches_, oI)
191  {
192  const label patchI = patches_[oI];
193 
194  tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
195  const vectorField& nf = tnf();
196  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
197 
198  bdJdvPtr_()[patchI] =
199  - (p.boundaryField()[patchI] + 0.5*magSqr(Ub))*nf
200  - (Ub & nf)*Ub;
201  }
202 }
203 
204 
206 {
207  const volScalarField& p = vars_.p();
208  const volVectorField& U = vars_.U();
209 
210  forAll(patches_, oI)
211  {
212  const label patchI = patches_[oI];
213 
214  tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
215  const vectorField& nf = tnf();
216 
217  bdJdvnPtr_()[patchI] =
218  - p.boundaryField()[patchI]
219  - 0.5*magSqr(U.boundaryField()[patchI])
220  - sqr(U.boundaryField()[patchI] & nf);
221  }
222 }
223 
224 
226 {
227  const volVectorField& U = vars_.U();
228 
229  forAll(patches_, oI)
230  {
231  const label patchI = patches_[oI];
232 
233  tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
234  const vectorField& nf = tnf();
235  scalarField Un(U.boundaryField()[patchI] & nf);
236 
237  bdJdvtPtr_()[patchI] = -Un*(U.boundaryField()[patchI] - Un*nf);
238  }
239 }
240 
241 
243 {
244  if (Pstream::master())
245  {
246  // file is opened only upon invocation of the write function
247  // in order to avoid various instantiations of the same objective
248  // opening the same file
249  unsigned int width = IOstream::defaultPrecision() + 5;
250  if (objFunctionFilePtr_.empty())
251  {
253  objFunctionFilePtr_() << setw(4) << "#" << " ";
254  objFunctionFilePtr_() << setw(width) << "ptLosses" << " ";
255  forAll(patches_, oI)
256  {
257  label patchI = patches_[oI];
259  << setw(width) << mesh_.boundary()[patchI].name() << " ";
260  }
262  }
263 
264  objFunctionFilePtr_() << setw(4) << mesh_.time().value() << " ";
265  objFunctionFilePtr_() << setw(width) << J_ << " ";
266  forAll(patchPt_, pI)
267  {
268  objFunctionFilePtr_() << setw(width) << patchPt_[pI] << " ";
269  }
271  }
272 }
273 
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 } // End namespace objectives
278 } // End namespace Foam
279 
280 // ************************************************************************* //
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::objectives::objectivePtLosses::write
void write() const
Write objective function values and its contrituents.
Definition: objectivePtLosses.C:242
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::DynamicList< label >
objectivePtLosses.H
Foam::objectiveIncompressible::vars_
const incompressibleVars & vars_
Definition: objectiveIncompressible.H:62
Foam::objectives::addToRunTimeSelectionTable
addToRunTimeSelectionTable(objectiveIncompressible, objectiveForce, dictionary)
Foam::wordReList
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:52
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::objective::setObjectiveFilePtr
void setObjectiveFilePtr() const
Set the output file ptr.
Definition: objective.C:58
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:185
Foam::objectiveIncompressible::bdJdvPtr_
autoPtr< boundaryVectorField > bdJdvPtr_
Definition: objectiveIncompressible.H:79
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
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
Foam::HashSet< label, Hash< label > >
Foam::objective::objFunctionFilePtr_
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:129
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
coupledFvPatch.H
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::objectives::defineTypeNameAndDebug
defineTypeNameAndDebug(objectiveForce, 0)
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< vector >
Foam::objective::J_
scalar J_
Definition: objective.H:73
createZeroField.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:436
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:121
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:84
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::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
U
U
Definition: pEqn.H:72
Foam::objectives::objectivePtLosses::J
scalar J()
Return the objective function value.
Definition: objectivePtLosses.C:140
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::objectives::objectivePtLosses::update_boundarydJdp
void update_boundarydJdp()
Update values to be added to the adjoint inlet velocity.
Definition: objectivePtLosses.C:169
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:325
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:64
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:857
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::objective::dict
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:93
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
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:205
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::objectives::objectivePtLosses::update_boundarydJdvt
void update_boundarydJdvt()
Update values to be added to the adjoint outlet tangential velocity.
Definition: objectivePtLosses.C:225
Foam::incompressibleVars::pInst
const volScalarField & pInst() const
Return const reference to pressure.
Definition: incompressibleVars.C:382