Loading...
Searching...
No Matches
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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
43(
45 boundaryAdjointContributionIncompressible,
47);
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52boundaryAdjointContributionIncompressible::
53boundaryAdjointContributionIncompressible
54(
55 const word& managerName,
57 const word& simulationType,
58 const fvPatch& patch
59)
60:
61 boundaryAdjointContribution
62 (
63 managerName,
65 simulationType,
66 patch
67 ),
69 (
70 patch_.patch().boundaryMesh().mesh().
71 lookupObjectRef<objectiveManager>(managerName)
72 ),
74 (
75 patch_.patch().boundaryMesh().mesh().
77 getPrimalVars()
78 ),
80 (
81 patch_.patch().boundaryMesh().mesh().
83 )
84{}
85
86
87// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88
90{
91 // Objective function contribution
92 tmp<vectorField> tsource =
94 (
98 );
99 vectorField& source = tsource.ref();
100
101 // Turbulence model differentiation contribution.
104 source += adjointRAS().adjointMomentumBCSource()[patch_.index()];
105
106 return tsource;
107}
108
109
111{
112 // Objective function contribution
113 tmp<scalarField> tsource =
115 (
119 );
120
121 scalarField& source = tsource.ref();
122
123 // Turbulence model differentiation contribution.
126 const vectorField& adjointTurbulenceContr =
127 adjointRAS().adjointMomentumBCSource()[patch_.index()];
128
129 source += adjointTurbulenceContr & patch_.nf();
130
131 return (tsource);
132}
133
134
137{
138 // Objective function contribution
139 tmp<vectorField> tsource =
141 (
145 );
146
147 vectorField& source = tsource.ref();
148
149 // Turbulence model differentiation contribution.
152 const vectorField& adjointTurbulenceContr =
153 adjointRAS().adjointMomentumBCSource()[patch_.index()];
154
155 tmp<vectorField> tnf = patch_.nf();
156 const vectorField& nf = tnf();
157
158 source += adjointTurbulenceContr - (adjointTurbulenceContr & nf)*nf;
159
160 return (tsource);
161}
162
163
166{
167 return
169 (
182 (
196 (
210 (
224 (
237 (
243
244
260 nu = turbulenceModel().nu()().boundaryField()[patch_.index()];
261
262 return tnu;
263}
264
265
267{
268 /*
269 const polyMesh& mesh = patch_.patch().boundaryMesh().mesh();
270 const compressible::turbulenceModel& turbulenceModel =
271 mesh.lookupObject<compressible::turbulenceModel>("turbulenceModel");
272 tmp<scalarField> talphaEff = turbulenceModel.alphaEff(patch_.index());
273 */
274
276 << "no abstract thermalDiffusion is implemented. Returning zero field";
277
278 return tmp<scalarField>::New(patch_.size(), Zero);
279}
280
281
283{
284 return primalVars_.turbulence()->y()[patch_.index()];
285}
286
287
290{
291 return
292 adjointVars().adjointTurbulence()->diffusionCoeffVar1(patch_.index());
293
294}
295
296
299{
300 return
301 adjointVars().adjointTurbulence()->diffusionCoeffVar2(patch_.index());
302}
303
304
307 return
308 primalVars_.RASModelVariables()->TMVar1().
309 boundaryField()[patch_.index()];
310}
311
312
315 return
316 primalVars_.RASModelVariables()->TMVar2().
317 boundaryField()[patch_.index()];
318}
319
322{
323 return primalVars_.U().boundaryField()[patch_.index()];
324}
325
326
328{
329 return primalVars_.p().boundaryField()[patch_.index()];
330}
331
332
335{
336 return primalVars_.phi().boundaryField()[patch_.index()];
337}
338
339
342{
343 return primalVars_.RASModelVariables()().nutPatchField(patch_.index());
344}
345
350}
351
352
357
358
363}
364
369}
370
371
376
377
383
384
387{
388 return adjointSolver_.getAdjointVars();
389}
390
391
394{
395 return objectiveManager_;
396}
397
398
399// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
400
401} // End namespace Foam
402
403// ************************************************************************* //
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 noexcept
Return const-reference to the boundary field.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Contributions of objective function differentiation to adjoint boundary conditions for incompressible...
const incompressibleAdjointSolver & adjointSolver_
Note: getting a reference to the adjoint vars in the constructor of boundaryAdjointContributionIncomp...
Abstract base class for computing contributions of the objective functions to the adjoint boundary co...
tmp< Field< returnType > > sumContributions(PtrList< sourceType > &sourceList, const fvPatchField< returnType > &(castType::*boundaryFunction)(const label), bool(castType::*hasFunction)() const)
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition fvPatch.H:242
label index() const noexcept
The index of this patch in the boundary mesh.
Definition fvPatch.H:218
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition fvPatch.C:143
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.
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 autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
const fvPatchScalarField & boundarydJdnut(const label)
Objective partial deriv wrt nut for a specific patch.
const fvPatchTensorField & boundarydJdGradU(const label)
Objective partial deriv wrt stress tensor.
const fvPatchVectorField & boundarydJdvt(const label)
Objective partial deriv wrt tangent velocity for a specific patch.
const fvPatchVectorField & boundarydJdv(const label)
Objective partial deriv wrt velocity for a specific patch.
const fvPatchScalarField & boundarydJdTMvar2(const label)
Objective partial deriv wrt turbulence model var 2 for a specific patch.
const fvPatchScalarField & boundarydJdvn(const label)
Objective partial deriv wrt normal velocity for a specific patch.
const fvPatchScalarField & boundarydJdTMvar1(const label)
Objective partial deriv wrt turbulence model var 1 for a specific patch.
const fvPatchScalarField & boundarydJdT(const label)
Objective partial deriv wrt temperature for a specific patch.
const fvPatchVectorField & boundarydJdp(const label)
Objective partial deriv wrt pressure (times normal) for a specific patch.
Class for managing objective functions.
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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.
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
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< vector > fvPatchVectorField
fvsPatchField< scalar > fvsPatchScalarField
fvPatchField< scalar > fvPatchScalarField
volScalarField & nu