objectiveForce.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 "objectiveForce.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38namespace objectives
39{
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
45(
49);
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const fvMesh& mesh,
57 const dictionary& dict,
58 const word& adjointSolverName,
59 const word& primalSolverName
60)
61:
62 objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
63 forcePatches_
64 (
65 mesh_.boundaryMesh().patchSet
66 (
67 dict.get<wordRes>("patches")
68 ).sortedToc()
69 ),
70 forceDirection_(dict.get<vector>("direction")),
71 Aref_(dict.get<scalar>("Aref")),
72 rhoInf_(dict.get<scalar>("rhoInf")),
73 UInf_(dict.get<scalar>("UInf"))
74{
75 // Sanity check and print info
76 if (forcePatches_.empty())
77 {
79 << "No valid patch name on which to minimize " << type() << endl
80 << exit(FatalError);
81 }
82 if (debug)
83 {
84 Info<< "Minimizing " << type() << " in patches:" << endl;
85 for (const label patchI : forcePatches_)
86 {
87 Info<< "\t " << mesh_.boundary()[patchI].name() << endl;
88 }
89 }
90
91 // Allocate boundary field pointers
92 bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
93 bdSdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
94 bdxdbMultPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
95 bdJdnutPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
96 bdJdGradUPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
97}
98
99
100// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101
103{
104 vector pressureForce(Zero);
105 vector viscousForce(Zero);
106 vector cumulativeForce(Zero);
107
108
109 const volScalarField& p = vars_.pInst();
112
113 volSymmTensorField devReff(turbulence->devReff());
114
115 for (const label patchI : forcePatches_)
116 {
117 const vectorField& Sf = mesh_.Sf().boundaryField()[patchI];
118 pressureForce += gSum(Sf*p.boundaryField()[patchI]);
119 viscousForce += gSum(devReff.boundaryField()[patchI] & Sf);
120 }
121
122 cumulativeForce = pressureForce + viscousForce;
123
124 scalar force = cumulativeForce & forceDirection_;
125
126 // Intentionally not using denom - derived might implement virtual denom()
127 // function differently
128 scalar Cforce = force/(0.5*UInf_*UInf_*Aref_);
129
131 << "Force|Coeff " << force << "|" << Cforce << endl;
132
133 J_ = Cforce;
134
135 return Cforce;
136}
137
138
140{
141 for (const label patchI : forcePatches_)
142 {
143 bdJdpPtr_()[patchI] = forceDirection_/denom();
144 }
145}
146
147
149{
150 // Compute contributions with mean fields, if present
151 const volScalarField& p = vars_.p();
152 const volVectorField& U = vars_.U();
156
157 tmp<volSymmTensorField> tdevReff = turbVars->devReff(lamTransp, U);
158 const volSymmTensorField& devReff = tdevReff();
159
160 for (const label patchI : forcePatches_)
161 {
162 bdSdbMultPtr_()[patchI] =
163 (
164 (
165 forceDirection_& devReff.boundaryField()[patchI]
166 )
167 + (forceDirection_)*p.boundaryField()[patchI]
168 )
169 /denom();
170 }
171}
172
173
175{
176 const volScalarField& p = vars_.p();
177 const volVectorField& U = vars_.U();
178
180 turbVars = vars_.RASModelVariables();
182
183 // We only need to modify the boundaryField of gradU locally.
184 // If grad(U) is cached then
185 // a. The .ref() call fails since the tmp is initialised from a
186 // const ref
187 // b. we would be changing grad(U) for all other places in the code
188 // that need it
189 // So, always allocate new memory and avoid registering the new field
190 tmp<volTensorField> tgradU =
191 volTensorField::New("gradULocal", fvc::grad(U));
192 volTensorField& gradU = tgradU.ref();
194
195 // Explicitly correct the boundary gradient to get rid of
196 // the tangential component
197 forAll(mesh_.boundary(), patchI)
198 {
199 const fvPatch& patch = mesh_.boundary()[patchI];
200 if (isA<wallFvPatch>(patch))
201 {
202 tmp<vectorField> tnf = patch.nf();
203 gradUbf[patchI] = tnf*U.boundaryField()[patchI].snGrad();
204 }
205 }
206
207 // Term coming from gradp
209 const volVectorField& gradp = tgradp.cref();
210 for (const label patchI : forcePatches_)
211 {
212 bdxdbMultPtr_()[patchI] =
213 (forceDirection_ & mesh_.boundary()[patchI].nf())
214 *gradp.boundaryField()[patchI]/denom();
215 }
216 tgradp.clear();
217
218 // Term coming from stresses
219 tmp<volScalarField> tnuEff = lamTransp.nu() + turbVars->nutRef();
220 tmp<volSymmTensorField> tstress = tnuEff*twoSymm(tgradU);
221 const volSymmTensorField& stress = tstress.cref();
223 (Foam::createZeroFieldPtr<vector>( mesh_, "temp", sqr(dimVelocity)));
224 volVectorField& temp = ptemp.ref();
225
226 for (label idir = 0; idir < pTraits<vector>::nComponents; ++idir)
227 {
228 unzipRow(stress, idir, temp);
229 volTensorField gradStressDir(fvc::grad(temp));
230 for (const label patchI : forcePatches_)
231 {
232 const fvPatch& patch = mesh_.boundary()[patchI];
233 tmp<vectorField> tnf = patch.nf();
234 bdxdbMultPtr_()[patchI] -=
235 forceDirection_.component(idir)
236 *(gradStressDir.boundaryField()[patchI] & tnf)/denom();
237 }
238 }
239}
240
241
243{
244 const volVectorField& U = vars_.U();
246
247 for (const label patchI : forcePatches_)
248 {
249 const fvPatch& patch = mesh_.boundary()[patchI];
250 tmp<vectorField> tnf = patch.nf();
251 bdJdnutPtr_()[patchI] =
252 - ((devGradU.boundaryField()[patchI] & forceDirection_) & tnf)
253 /denom();
254 }
255}
256
257
259{
263 volScalarField nuEff(lamTransp.nu() + turbVars->nutRef());
264 for (const label patchI : forcePatches_)
265 {
266 const fvPatch& patch = mesh_.boundary()[patchI];
267 const vectorField& Sf = patch.Sf();
268 bdJdGradUPtr_()[patchI] =
269 - nuEff.boundaryField()[patchI]
271 }
272}
273
274
276{
277 return 0.5*UInf_*UInf_*Aref_;
278}
279
280
282{
283 return forceDirection_;
284}
285
286
287// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288
289} // End namespace objectives
290} // End namespace Foam
291
292// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
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
const surfaceVectorField & Sf() const
Return cell face area vectors.
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< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
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
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
const vector & forceDirection() const
Return force direction.
void update_boundarydJdnut()
Update dJ/dnut multiplier.
virtual scalar denom() const
Return denominator, without density.
void update_boundarydJdGradU()
Update dJ/dGradU 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