boundaryAdjointContributionIncompressible.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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 FOSS GP
10 Copyright (C) 2019 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
31#include "adjointRASModel.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
43(
47);
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52boundaryAdjointContributionIncompressible::
53boundaryAdjointContributionIncompressible
54(
55 const word& managerName,
56 const word& adjointSolverName,
57 const word& simulationType,
58 const fvPatch& patch
59)
60:
62 (
63 managerName,
64 adjointSolverName,
65 simulationType,
66 patch
67 ),
68 objectiveManager_
69 (
70 patch_.patch().boundaryMesh().mesh().
71 lookupObjectRef<objectiveManager>(managerName)
72 ),
73 primalVars_
74 (
75 patch_.patch().boundaryMesh().mesh().
76 lookupObject<incompressibleAdjointSolver>(adjointSolverName).
77 getPrimalVars()
78 ),
79 adjointSolver_
80 (
81 patch_.patch().boundaryMesh().mesh().
82 lookupObject<incompressibleAdjointSolver>(adjointSolverName)
83 )
84{}
85
86
87// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88
90{
91 // Objective function contribution
92 tmp<vectorField> tsource =
94 (
97 );
98 vectorField& source = tsource.ref();
99
100 // Turbulence model differentiation contribution.
103 source += adjointRAS().adjointMomentumBCSource()[patch_.index()];
104
105 return tsource;
106}
107
108
110{
111 // Objective function contribution
112 tmp<scalarField> tsource =
114 (
117 );
118
119 scalarField& source = tsource.ref();
120
121 // Turbulence model differentiation contribution.
124 const vectorField& adjointTurbulenceContr =
125 adjointRAS().adjointMomentumBCSource()[patch_.index()];
126
127 source += adjointTurbulenceContr & patch_.nf();
128
129 return (tsource);
130}
131
132
135{
136 // Objective function contribution
137 tmp<vectorField> tsource =
139 (
142 );
143
144 vectorField& source = tsource.ref();
145
146 // Turbulence model differentiation contribution.
149 const vectorField& adjointTurbulenceContr =
150 adjointRAS().adjointMomentumBCSource()[patch_.index()];
151
152 tmp<vectorField> tnf = patch_.nf();
153 const vectorField& nf = tnf();
154
155 source += adjointTurbulenceContr - (adjointTurbulenceContr & nf)*nf;
156
157 return (tsource);
158}
159
160
163{
164 return
166 (
169 );
170
171}
172
173
175{
176 return
178 (
181 );
182
183}
184
185
188{
189 return
191 (
194 );
195
196}
197
198
201{
202 return
204 (
207 );
208
209}
210
211
214{
215 return
217 (
220 );
221}
222
223
226{
227 return
229 (
232 );
233}
234
235
237{
238
239 return adjointVars().adjointTurbulence()().nuEff(patch_.index());
240}
241
242
244{
246 scalarField& nu = tnu.ref();
247
250
251 nu = turbulenceModel().nu()().boundaryField()[patch_.index()];
252
253 return tnu;
254}
255
256
258{
259 /*
260 const polyMesh& mesh = patch_.patch().boundaryMesh().mesh();
261 const compressible::turbulenceModel& turbulenceModel =
262 mesh.lookupObject<compressible::turbulenceModel>("turbulenceModel");
263 tmp<scalarField> talphaEff = turbulenceModel.alphaEff(patch_.index());
264 */
265
266 tmp<scalarField> talphaEff(new scalarField(patch_.size(), Zero));
267
269 << "no abstract thermalDiffusion is implemented. Returning zero field";
270
271
272 return talphaEff;
273}
274
275
277{
278 return primalVars_.turbulence()->y()[patch_.index()];
279}
280
281
284{
285 return
286 adjointVars().adjointTurbulence()->diffusionCoeffVar1(patch_.index());
287
288}
289
290
293{
294 return
295 adjointVars().adjointTurbulence()->diffusionCoeffVar2(patch_.index());
296}
297
298
300{
301 return
302 primalVars_.RASModelVariables()->TMVar1().
303 boundaryField()[patch_.index()];
304}
305
306
308{
309 return
310 primalVars_.RASModelVariables()->TMVar2().
311 boundaryField()[patch_.index()];
312}
313
314
316{
317 return primalVars_.U().boundaryField()[patch_.index()];
318}
319
320
322{
323 return primalVars_.p().boundaryField()[patch_.index()];
324}
325
326
329{
331}
332
333
336{
337 return
338 primalVars_.RASModelVariables()().nutRef().boundaryField()
339 [
340 patch_.index()
341 ];
342}
343
344
346{
348}
349
350
352{
354}
355
356
359{
361}
362
363
365{
366 return primalVars_.solverName();
367}
368
369
371{
372 return adjointVars().solverName();
373}
374
375
378{
379 return primalVars_;
380}
381
382
385{
387}
388
389
392{
393 return objectiveManager_;
394}
395
396
397// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398
399} // End namespace Foam
400
401// ************************************************************************* //
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.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Contributions of objective function differentiation to adjoint boundary conditions for incompressible...
tmp< Field< returnType > > sumContributions(PtrList< sourceType > &sourceList, const fvPatchField< returnType > &(castType::*boundaryFunction)(const label))
Abstract base class for computing contributions of the objective functions to the adjoint boundary co...
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual label size() const
Return size.
Definition: fvPatch.H:185
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:203
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:144
const volVectorField & UaInst() const
Return const reference to velocity.
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
const volScalarField & paInst() const
Return const reference to pressure.
Base class for incompressibleAdjoint solvers.
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
Class including all adjoint fields for incompressible flows.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
Base class for solution control classes.
const volVectorField & U() const
Return const reference to velocity.
const volScalarField & p() const
Return const reference to pressure.
const surfaceScalarField & phi() const
Return const reference to volume flux.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
const boundaryVectorField & boundarydJdp()
Objective partial deriv wrt pressure (times normal) for all patches.
const boundaryScalarField & boundarydJdTMvar1()
Objective partial deriv wrt turbulence model var 1 for all patches.
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
const boundaryTensorField & boundarydJdGradU()
Objective partial deriv wrt gradU.
const boundaryScalarField & boundarydJdnut()
Objective partial deriv wrt nut for all patches.
const boundaryVectorField & boundarydJdv()
Objective partial deriv wrt velocity for all patches.
const boundaryScalarField & boundarydJdvn()
Objective partial deriv wrt normal velocity for all patches.
const boundaryVectorField & boundarydJdvt()
Objective partial deriv wrt tangent velocity for all patches.
const boundaryScalarField & boundarydJdTMvar2()
Objective partial deriv wrt turbulence model var 2 for all patches.
class for managing incompressible objective functions.
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
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
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
volScalarField & nu