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 -------------------------------------------------------------------------------
12 License
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 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(boundaryAdjointContributionIncompressible, 0);
43 (
44  boundaryAdjointContribution,
45  boundaryAdjointContributionIncompressible,
46  dictionary
47 );
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 boundaryAdjointContributionIncompressible::
53 boundaryAdjointContributionIncompressible
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  (
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  (
173  );
174 
175 }
176 
177 
179 {
180  return
182  (
186  );
187 
188 }
189 
190 
193 {
194  return
196  (
200  );
201 
202 }
203 
204 
207 {
208  return
210  (
214  );
215 
216 }
217 
218 
221 {
222  return
224  (
228  );
229 }
230 
231 
234 {
235  return
237  (
241  );
242 }
243 
244 
246 {
247 
248  return adjointVars().adjointTurbulence()().nuEff(patch_.index());
249 }
250 
251 
253 {
255  scalarField& nu = tnu.ref();
256 
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 
275  tmp<scalarField> talphaEff(new scalarField(patch_.size(), Zero));
276 
278  << "no abstract thermalDiffusion is implemented. Returning zero field";
279 
280 
281  return talphaEff;
282 }
283 
284 
286 {
287  return primalVars_.turbulence()->y()[patch_.index()];
288 }
289 
290 
293 {
294  return
295  adjointVars().adjointTurbulence()->diffusionCoeffVar1(patch_.index());
296 
297 }
298 
299 
302 {
303  return
304  adjointVars().adjointTurbulence()->diffusionCoeffVar2(patch_.index());
305 }
306 
307 
309 {
310  return
311  primalVars_.RASModelVariables()->TMVar1().
312  boundaryField()[patch_.index()];
313 }
314 
315 
317 {
318  return
319  primalVars_.RASModelVariables()->TMVar2().
320  boundaryField()[patch_.index()];
321 }
322 
325 {
326  return primalVars_.U().boundaryField()[patch_.index()];
327 }
328 
329 
331 {
332  return primalVars_.p().boundaryField()[patch_.index()];
333 }
334 
335 
336 const fvsPatchScalarField&
338 {
339  return primalVars_.phi().boundaryField()[patch_.index()];
340 }
341 
342 
345 {
346  return primalVars_.RASModelVariables()().nutPatchField(patch_.index());
347 }
348 
351 {
352  return adjointVars().UaInst().boundaryField()[patch_.index()];
353 }
354 
355 
357 {
358  return adjointVars().paInst().boundaryField()[patch_.index()];
359 }
360 
361 
364 {
366 }
367 
370 {
371  return primalVars_.solverName();
372 }
373 
374 
376 {
377  return adjointVars().solverName();
378 }
379 
380 
381 const incompressibleVars&
383 {
384  return primalVars_;
385 }
386 
387 
390 {
392 }
393 
394 
397 {
398  return objectiveManager_;
399 }
400 
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 } // End namespace Foam
405 
406 // ************************************************************************* //
const volScalarField & paInst() const
Return const reference to pressure.
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
bool hasBoundarydJdTMVar1() const noexcept
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition: fvPatch.C:143
Class for managing objective functions.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
Class including all adjoint fields for incompressible flows.
Base class for incompressibleAdjoint solvers.
const boundaryScalarField & boundarydJdTMvar1()
Objective partial deriv wrt turbulence model var 1 for all patches.
const surfaceScalarField & phi() const
Return const reference to volume flux.
Base class for solution control classes.
const boundaryVectorField & boundarydJdv()
Objective partial deriv wrt velocity for all patches.
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const boundaryScalarField & boundarydJdTMvar2()
Objective partial deriv wrt turbulence model var 2 for all patches.
const volScalarField & p() const
Return const reference to pressure.
const incompressibleAdjointSolver & adjointSolver_
Note: getting a reference to the adjoint vars in the constructor of boundaryAdjointContributionIncomp...
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
const word & solverName() const
Return solver name.
Definition: variablesSet.C:77
const volVectorField & U() const
Return const reference to velocity.
Abstract base class for computing contributions of the objective functions to the adjoint boundary co...
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition: fvPatch.H:234
defineTypeNameAndDebug(combustionModel, 0)
const boundaryTensorField & boundarydJdGradU()
Objective partial deriv wrt gradU.
bool hasBoundarydJdTMVar2() const noexcept
const boundaryScalarField & boundarydJdnut()
Objective partial deriv wrt nut for all patches.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
const volVectorField & UaInst() const
Return const reference to velocity.
const boundaryVectorField & boundarydJdvt()
Objective partial deriv wrt tangent velocity for all patches.
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
#define WarningInFunction
Report a warning using Foam::Warning.
const std::string patch
OpenFOAM patch number as a std::string.
const boundaryVectorField & boundarydJdp()
Objective partial deriv wrt pressure (times normal) for all patches.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
const boundaryScalarField & boundarydJdvn()
Objective partial deriv wrt normal velocity for all patches.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< Field< returnType > > sumContributions(PtrList< sourceType > &sourceList, const fvPatchField< returnType > &(castType::*boundaryFunction)(const label), bool(castType::*hasFunction)() const)
label index() const noexcept
The index of this patch in the boundary mesh.
Definition: fvPatch.H:218
volScalarField & nu
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127