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,
50  dictionary
51 );
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void objectiveUniformityPatch::initialize()
57 {
58  // If patches are prescribed, use them
59 
60  wordRes patchSelection;
61  if (dict().readIfPresent("patches", patchSelection))
62  {
63  patches_ = mesh_.boundaryMesh().patchSet(patchSelection).sortedToc();
64  }
65  // Otherwise, pick them up based on the mass flow.
66  // Note: a non-zero U initialisation should be used in order to pick up the
67  // outlet patches correctly
68  else
69  {
71  << "No patches provided to " << type() << ". "
72  << "Choosing them according to the patch mass flows" << nl;
73 
74  DynamicList<label> objectiveReportPatches(mesh_.boundary().size());
76  forAll(mesh_.boundary(), patchI)
77  {
78  const fvsPatchScalarField& phiPatch = phi.boundaryField()[patchI];
79  if (!isA<coupledFvPatch>(mesh_.boundary()[patchI]))
80  {
81  const scalar mass = gSum(phiPatch);
82  if (mass > SMALL)
83  {
84  objectiveReportPatches.push_back(patchI);
85  }
86  }
87  }
88  patches_.transfer(objectiveReportPatches);
89  }
90  UMean_.setSize(patches_.size(), Zero);
91  UVar_.setSize(patches_.size(), Zero);
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;
105  }
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_(),
122  UMean_(),
123  UVar_()
124 {
125  // Find inlet/outlet patches
126  initialize();
127 
128  // Allocate boundary field pointers
129  bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
130  bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
131  bdJdvtPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
132 }
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  J_ = Zero;
140 
141  const volVectorField& U = vars_.UInst();
142 
143  forAll(patches_, oI)
144  {
145  const label patchI = patches_[oI];
146  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
147  const scalar sumMagSf = gSum(magSf);
148  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
149  UMean_[oI] = gSum(Ub*magSf)/sumMagSf;
150  UVar_[oI] = gSum(magSqr(Ub - UMean_[oI])*magSf)/sumMagSf;
151  J_ += 0.5*UVar_[oI];
152  }
153 
154  return J_;
155 }
156 
157 
159 {
160  const volVectorField& U = vars_.U();
161 
162  forAll(patches_, oI)
163  {
164  const label patchI = patches_[oI];
165  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
166  const scalar sumMagSf = gSum(magSf);
167  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
168  bdJdvPtr_()[patchI] = (Ub - UMean_[oI])/sumMagSf;
169  }
170 }
171 
172 
174 {
175  const volVectorField& U = vars_.U();
176 
177  forAll(patches_, oI)
178  {
179  const label patchI = patches_[oI];
180  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
181  const scalar sumMagSf = gSum(magSf);
182  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
183  tmp<vectorField> nf = mesh_.boundary()[patchI].nf();
184  bdJdvnPtr_()[patchI] = ((Ub - UMean_[oI]) & nf)/sumMagSf;
185  }
186 }
187 
188 
190 {
191  const volVectorField& U = vars_.U();
192 
193  forAll(patches_, oI)
194  {
195  const label patchI = patches_[oI];
196  const scalarField& magSf = mesh_.boundary()[patchI].magSf();
197  const scalar sumMagSf = gSum(magSf);
198  const fvPatchVectorField& Ub = U.boundaryField()[patchI];
199  tmp<vectorField> nf = mesh_.boundary()[patchI].nf();
200  vectorField UdiffTangent(Ub - UMean_[oI]);
202  bdJdvtPtr_()[patchI] =
203  (UdiffTangent - (UdiffTangent & nf())*nf())/sumMagSf;
204  }
205 }
206 
207 
209 {
210  OFstream& file = objFunctionFilePtr_();
211  for (const label patchI : patches_)
212  {
213  const word patchName(mesh_.boundary()[patchI].name());
214  file<< setw(width_) << word(patchName + "-" + "UMean") << " ";
215  file<< setw(width_) << word(patchName + "-" + "UVar") << " ";
216  file<< setw(width_) << word(patchName + "-" + "UStd") << " ";
217  }
218 }
219 
220 
222 {
223  OFstream& file = objFunctionFilePtr_();
224  forAll(UMean_, pI)
225  {
226  file<< setw(width_) << mag(UMean_[pI]) << " ";
227  file<< setw(width_) << UVar_[pI] << " ";
228  file<< setw(width_) << sqrt(UVar_[pI]) << " ";
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 } // End namespace objectives
236 } // End namespace Foam
237 
238 // ************************************************************************* //
fvsPatchField< scalar > fvsPatchScalarField
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:91
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:326
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:209
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.
Output to file stream, using an OSstream.
Definition: OFstream.H:49
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
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:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
defineTypeNameAndDebug(objectivePartialVolume, 1)
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const fvMesh & mesh_
Definition: objective.H:63
scalar J_
Objective function value and weight.
Definition: objective.H:76
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...
addToRunTimeSelectionTable(objectiveGeometric, objectivePartialVolume, dictionary)
virtual 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:226
virtual void update_boundarydJdv()
Update values to be added to the adjoint outlet velocity.
virtual 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:78
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual void update_boundarydJdvn()
Update values to be added to the adjoint outlet pressure.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const volVectorField & UInst() const
Return const reference to velocity.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
Abstract base class for objective functions in incompressible flows.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127