objectiveUniformityPatch.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 -------------------------------------------------------------------------------
11 License
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 
30 #include "createZeroField.H"
31 #include "coupledFvPatch.H"
32 #include "IOmanip.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 namespace objectives
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(objectiveUniformityPatch, 0);
47 (
48  objectiveIncompressible,
49  objectiveUniformityPatch,
51 );
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const fvMesh& mesh,
59  const dictionary& dict,
60  const word& adjointSolverName,
61  const word& primalSolverName
62 )
63 :
64  objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
65  patches_(),
66  UMean_(),
67  UVar_()
68 {
69  // Find inlet/outlet patches
70  initialize();
71 
72  // Allocate boundary field pointers
73  bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
74  bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
75  bdJdvtPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
76 }
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82 {
83  // If patches are prescribed, use them
84 
85  wordRes patchSelection;
86  if (dict().readIfPresent("patches", patchSelection))
87  {
88  patches_ = mesh_.boundaryMesh().patchSet(patchSelection).sortedToc();
89  }
90  // Otherwise, pick them up based on the mass flow.
91  // Note: a non-zero U initialisation should be used in order to pick up the
92  // outlet patches correctly
93  else
94  {
96  << "No patches provided to " << type() << ". "
97  << "Choosing them according to the patch mass flows" << nl;
98 
99  DynamicList<label> objectiveReportPatches(mesh_.boundary().size());
100  const surfaceScalarField& phi = vars_.phiInst();
101  forAll(mesh_.boundary(), patchI)
102  {
103  const fvsPatchScalarField& phiPatch = phi.boundaryField()[patchI];
104  if (!isA<coupledFvPatch>(mesh_.boundary()[patchI]))
105  {
106  const scalar mass = gSum(phiPatch);
107  if (mass > SMALL)
108  {
109  objectiveReportPatches.append(patchI);
110  }
111  }
112  }
113  patches_.transfer(objectiveReportPatches);
114  }
115  UMean_.setSize(patches_.size(), Zero);
116  UVar_.setSize(patches_.size(), Zero);
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  const volVectorField& U = vars_.UInst();
140 
141  forAll(patches_, oI)
142  {
143  const label patchI = patches_[oI];
144  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
145  const scalar sumMagSf = gSum(magSf);
146  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
147  UMean_[oI] = gSum(Ub*magSf)/sumMagSf;
148  UVar_[oI] = gSum(magSqr(Ub - UMean_[oI])*magSf)/sumMagSf;
149  J_ += 0.5*UVar_[oI];
150  }
151 
152  return J_;
153 }
154 
155 
157 {
158  const volVectorField& U = vars_.U();
159 
160  forAll(patches_, oI)
161  {
162  const label patchI = patches_[oI];
163  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
164  const scalar sumMagSf = gSum(magSf);
165  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
166  bdJdvPtr_()[patchI] = (Ub - UMean_[oI])/sumMagSf;
167  }
168 }
169 
170 
172 {
173  const volVectorField& U = vars_.U();
174 
175  forAll(patches_, oI)
176  {
177  const label patchI = patches_[oI];
178  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
179  const scalar sumMagSf = gSum(magSf);
180  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
181  tmp<vectorField> nf = mesh_.boundary()[patchI].nf();
182  bdJdvnPtr_()[patchI] = ((Ub - UMean_[oI]) & nf)/sumMagSf;
183  }
184 }
185 
186 
188 {
189  const volVectorField& U = vars_.U();
190 
191  forAll(patches_, oI)
192  {
193  const label patchI = patches_[oI];
194  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
195  const scalar sumMagSf = gSum(magSf);
196  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
197  tmp<vectorField> nf = mesh_.boundary()[patchI].nf();
198  vectorField UdiffTangent(Ub - UMean_[oI]);
200  bdJdvtPtr_()[patchI] =
201  (UdiffTangent - (UdiffTangent & nf())*nf())/sumMagSf;
202  }
203 }
204 
205 
207 {
208  OFstream& file = objFunctionFilePtr_();
209  for (const label patchI : patches_)
210  {
211  const word patchName(mesh_.boundary()[patchI].name());
212  file<< setw(width_) << word(patchName + "-" + "UMean") << " ";
213  file<< setw(width_) << word(patchName + "-" + "UVar") << " ";
214  file<< setw(width_) << word(patchName + "-" + "UStd") << " ";
215  }
216 }
217 
218 
220 {
221  OFstream& file = objFunctionFilePtr_();
222  forAll(UMean_, pI)
223  {
224  file<< setw(width_) << mag(UMean_[pI]) << " ";
225  file<< setw(width_) << UVar_[pI] << " ";
226  file<< setw(width_) << sqrt(UVar_[pI]) << " ";
227  }
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace objectives
234 } // End namespace Foam
235 
236 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
fvPatchField< vector > fvPatchVectorField
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:87
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:439
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:173
Output to file stream, using an OSstream.
Definition: OFstream.H:49
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:420
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
defineTypeNameAndDebug(objectiveFlowRate, 0)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
virtual void addHeaderColumns() const
Write headers for additional columns.
Macros for easy insertion into run-time selection tables.
objectiveUniformityPatch(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
from components
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
void initialize()
Return the objectiveReportPatches.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
const fvMesh & mesh_
Definition: objective.H:63
scalar J_
Objective function value and weight.
Definition: objective.H:75
autoPtr< boundaryVectorField > bdJdvPtr_
const volVectorField & U() const
Return const reference to velocity.
Istream and Ostream manipulators taking arguments.
int debug
Static debugging option.
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...
void update_boundarydJdvt()
Update values to be added to the adjoint outlet tangential velocity.
const incompressibleVars & vars_
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:190
void update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
addToRunTimeSelectionTable(objectiveIncompressible, objectiveFlowRate, dictionary)
scalar J()
Return the objective function value.
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
virtual void addColumnValues() const
Write information to additional columns.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:130
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
messageStream Info
Information stream (stdout output on master, null elsewhere)
void update_boundarydJdvn()
Update values to be added to the adjoint outlet pressure.
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.
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const volVectorField & UInst() const
Return const reference to velocity.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
Abstract base class for objective functions in incompressible flows.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:705
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157