objectiveMoment.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, 2022 PCOpt/NTUA
9 Copyright (C) 2013-2019, 2022 FOSS GP
10 Copyright (C) 2019-2020 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
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 "objectiveMoment.H"
31#include "createZeroField.H"
32#include "wallFvPatch.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40namespace objectives
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
47(
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 momentPatches_
66 (
67 mesh_.boundaryMesh().patchSet
68 (
69 dict.get<wordRes>("patches")
70 ).sortedToc()
71 ),
72 momentDirection_(dict.get<vector>("direction")),
73 rotationCentre_(dict.get<vector>("rotationCenter")),
74 Aref_(dict.get<scalar>("Aref")),
75 lRef_(dict.get<scalar>("lRef")),
76 rhoInf_(dict.get<scalar>("rhoInf")),
77 UInf_(dict.get<scalar>("UInf")),
78 invDenom_(2./(rhoInf_*UInf_*UInf_*Aref_*lRef_)),
79 devReff_(vars_.turbulence()->devReff()())
80{
81 // Sanity check and print info
82 if (momentPatches_.empty())
83 {
85 << "No valid patch name on which to minimize " << type() << endl
86 << exit(FatalError);
87 }
88 if (debug)
89 {
90 Info<< "Minimizing " << type() << " in patches:" << endl;
91 for (const label patchI : momentPatches_)
92 {
93 Info<< "\t " << mesh_.boundary()[patchI].name() << endl;
94 }
95 }
96
97 // Allocate boundary field pointers
98 bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
99 bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
100 bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
101 bdxdbDirectMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
102 bdJdnutPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
103 //bdJdGradUPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
104}
105
106
107// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108
110{
111 vector pressureMoment(Zero);
112 vector viscousMoment(Zero);
113 vector cumulativeMoment(Zero);
114
115 // Update field here and use the same value for all functions
116 const volScalarField& p = vars_.pInst();
117 devReff_ = vars_.turbulence()->devReff()();
118
119 for (const label patchI : momentPatches_)
120 {
121 const fvPatch& patch = mesh_.boundary()[patchI];
122 const vectorField& Sf = patch.Sf();
123 vectorField dx(patch.Cf() - rotationCentre_);
124 pressureMoment += gSum
125 (
126 rhoInf_*(dx ^ Sf)*p.boundaryField()[patchI]
127 );
128
129 // Viscous term calculated using the full tensor derivative
130 viscousMoment += gSum
131 (
132 rhoInf_*(dx^(devReff_.boundaryField()[patchI] & Sf))
133 );
134 }
135
136 cumulativeMoment = pressureMoment + viscousMoment;
137
138 scalar moment = cumulativeMoment & momentDirection_;
139 scalar Cm = moment*invDenom_;
140 DebugInfo<<
141 "Moment|Coeff " << moment << "|" << Cm << endl;
142 J_ = Cm;
143 return Cm;
144}
145
146
148{
150 {
151 const volVectorField& U = vars_.U();
155
156 devReff_ = turbVars->devReff(lamTransp, U)();
157 }
158}
159
160
162{
163 for (const label patchI : momentPatches_)
164 {
165 const fvPatch& patch = mesh_.boundary()[patchI];
166 tmp<vectorField> tdx = patch.Cf() - rotationCentre_;
167 bdJdpPtr_()[patchI] = (momentDirection_ ^ tdx)*invDenom_*rhoInf_;
168 }
169}
170
171
173{
174 const volScalarField& p = vars_.p();
175
176 for (const label patchI : momentPatches_)
177 {
178 const fvPatch& patch = mesh_.boundary()[patchI];
179 tmp<vectorField> tdx = patch.Cf() - rotationCentre_;
180 bdSdbMultPtr_()[patchI] =
181 (
182 (
183 rhoInf_*
184 (
185 (momentDirection_ ^ tdx()) &
186 (
187 devReff_.boundaryField()[patchI]
188 )
189 )
190 )
191 + rhoInf_*(momentDirection_ ^ tdx())*p.boundaryField()[patchI]
192 )
193 *invDenom_;
194 }
195}
196
197
199{
200 const volScalarField& p = vars_.p();
201 const volVectorField& U = vars_.U();
202
206
207 // We only need to modify the boundaryField of gradU locally.
208 // If grad(U) is cached then
209 // a. The .ref() call fails since the tmp is initialised from a
210 // const ref
211 // b. we would be changing grad(U) for all other places in the code
212 // that need it
213 // So, always allocate new memory and avoid registering the new field
214 tmp<volTensorField> tgradU =
215 volTensorField::New("gradULocal", fvc::grad(U));
216 volTensorField::Boundary& gradUbf = tgradU.ref().boundaryFieldRef();
217
218 // Explicitly correct the boundary gradient to get rid of the
219 // tangential component
220 forAll(mesh_.boundary(), patchI)
221 {
222 const fvPatch& patch = mesh_.boundary()[patchI];
223 if (isA<wallFvPatch>(patch))
224 {
225 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
226 gradUbf[patchI] = tnf*U.boundaryField()[patchI].snGrad();
227 }
228 }
229
230 // Term coming from gradp
232 const volVectorField& gradp = tgradp.cref();
233 for (const label patchI : momentPatches_)
234 {
235 const fvPatch& patch = mesh_.boundary()[patchI];
236 tmp<vectorField> tnf = patch.nf();
237 tmp<vectorField> tdx = patch.Cf() - rotationCentre_;
238 bdxdbMultPtr_()[patchI] =
239 (momentDirection_ & (tdx ^ tnf))*gradp.boundaryField()[patchI]
240 *invDenom_*rhoInf_;
241 }
242 tgradp.clear();
243
244 // Term coming from stresses
245 tmp<volScalarField> tnuEff = lamTransp.nu() + turbVars->nutRef();
246 tmp<volSymmTensorField> tstress = tnuEff*twoSymm(tgradU);
247 const volSymmTensorField& stress = tstress.cref();
249 (Foam::createZeroFieldPtr<vector>( mesh_, "temp", sqr(dimVelocity)));
250 volVectorField& temp = ptemp.ref();
251
252 for (label idir = 0; idir < pTraits<vector>::nComponents; ++idir)
253 {
254 unzipRow(stress, idir, temp);
255 volTensorField gradStressDir(fvc::grad(temp));
256 for (const label patchI : momentPatches_)
257 {
258 const fvPatch& patch = mesh_.boundary()[patchI];
259 tmp<vectorField> tnf = patch.nf();
260 tmp<vectorField> tdx = patch.Cf() - rotationCentre_;
261 tmp<scalarField> taux = (momentDirection_ ^ tdx)().component(idir);
262 bdxdbMultPtr_()[patchI] -=
263 taux*(gradStressDir.boundaryField()[patchI] & tnf)
264 *invDenom_*rhoInf_;
265 }
266 }
267}
268
269
271{
272 const volScalarField& p = vars_.p();
273
274 for (const label patchI : momentPatches_)
275 {
276 const fvPatch& patch = mesh_.boundary()[patchI];
277 tmp<vectorField> tnf = patch.nf();
278 const vectorField& nf = tnf();
279 const vectorField dx(patch.Cf() - rotationCentre_);
280 const vectorField force
281 (
282 rhoInf_
283 *(
284 ((p.boundaryField()[patchI]*nf)
285 + (devReff_.boundaryField()[patchI] & nf))
286 )
287 );
288 bdxdbDirectMultPtr_()[patchI] =
289 (force^momentDirection_)*invDenom_*rhoInf_;
290 }
291}
292
293
295{
296 const volVectorField& U = vars_.U();
298
299 for (const label patchI : momentPatches_)
300 {
301 const fvPatch& patch = mesh_.boundary()[patchI];
302 tmp<vectorField> tnf = patch.nf();
303 tmp<vectorField> tdx = patch.Cf() - rotationCentre_;
304 const fvPatchSymmTensorField& bdevGradU =
305 devGradU.boundaryField()[patchI];
306 bdJdnutPtr_()[patchI] =
307 - rhoInf_*((tdx ^ (bdevGradU & tnf)) & momentDirection_)*invDenom_;
308 }
309}
310
311
312// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313
314} // End namespace objectives
315} // End namespace Foam
316
317// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTableI.H:59
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
T & ref()
Return reference to the managed object without nullptr checking.
Definition: autoPtr.H:161
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const volVectorField & U() const
Return const reference to velocity.
const volScalarField & p() const
Return const reference to pressure.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
const volScalarField & pInst() const
Return const reference to pressure.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Abstract base class for objective functions in incompressible flows.
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
const incompressibleVars & vars_
const fvMesh & mesh_
Definition: objective.H:67
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:109
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Definition: objective.H:120
bool computeMeanFields_
Definition: objective.H:72
scalar J_
Objective function value and weight.
Definition: objective.H:77
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:115
void update_meanValues()
Update mean drag and lift values.
void update_boundarydJdnut()
Update dJ/dnut multiplier.
scalar J()
Return the objective function value.
void update_dxdbMultiplier()
Update delta(x)/delta b multiplier.
void update_boundarydJdp()
Update values to be added to the adjoint wall velocity.
void update_dSdbMultiplier()
Update delta(n dS)/delta b multiplier.
A simple single-phase transport model based on viscosityModel.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
A class for managing temporary objects.
Definition: tmp.H:65
const T & cref() const
Definition: tmpI.H:213
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
compressible::turbulenceModel & turbulence
#define DebugInfo
Report an information message using Foam::Info.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
void unzipRow(const FieldField< Field, SymmTensor< Cmpt > > &input, const direction idx, FieldField< Field, Vector< Cmpt > > &result)
Extract a symmTensor field field row (x,y,z) == (0,1,2)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333