Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
37namespace Foam
40namespace objectives
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
47(
51);
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
56void 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;
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
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{
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{
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// ************************************************************************* //
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
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
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 surfaceScalarField & phiInst() const
Return const reference to volume flux.
Abstract base class for objective functions in incompressible flows.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
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
Computes and minimizes (half) the variance of the velocity distribution in a given set of patches.
objectiveUniformityPatch(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_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
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
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
dimensionedScalar sqrt(const dimensionedScalar &ds)
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