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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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
93  tmp<vectorField> tsource =
95  (
96  objectives,
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
114  tmp<scalarField> tsource =
116  (
117  objectives,
119  );
120 
121  scalarField& source = tsource.ref();
122 
123  // Turbulence model differentiation contribution.
126  const vectorField& adjointTurbulenceContr =
127  adjointRAS().adjointMomentumBCSource()[patch_.index()];
128 
129  tmp<vectorField> tnf = patch_.nf();
130  const vectorField& nf = tnf();
131 
132  source += adjointTurbulenceContr & nf;
133 
134  return (tsource);
135 }
136 
137 
140 {
141  // Objective function contribution
143  tmp<vectorField> tsource =
145  (
146  objectives,
148  );
149 
150  vectorField& source = tsource.ref();
151 
152  // Turbulence model differentiation contribution.
155  const vectorField& adjointTurbulenceContr =
156  adjointRAS().adjointMomentumBCSource()[patch_.index()];
157 
158  tmp<vectorField> tnf = patch_.nf();
159  const vectorField& nf = tnf();
160 
161  source += adjointTurbulenceContr - (adjointTurbulenceContr & nf)*nf;
162 
163  return (tsource);
164 }
165 
166 
169 {
171  tmp<vectorField> tsource =
173  (
174  objectives,
176  );
177 
178  return (tsource);
179 }
180 
181 
183 {
185  tmp<scalarField> tsource =
187  (
188  objectives,
190  );
191 
192  return (tsource);
193 }
194 
195 
198 {
200  tmp<scalarField> tsource =
202  (
203  objectives,
205  );
206 
207  return (tsource);
208 }
209 
210 
213 {
215  tmp<scalarField> tsource =
217  (
218  objectives,
220  );
221 
222  return (tsource);
223 }
224 
225 
227 {
228  tmp<scalarField> tnuEff(new scalarField(patch_.size(), Zero));
229  scalarField& nuEff = tnuEff.ref();
230 
232  adjointTurbulenceModel = adjointVars().adjointTurbulence();
233 
234  nuEff = adjointTurbulenceModel().nuEff()().boundaryField()[patch_.index()];
235 
236  return tnuEff;
237 }
238 
239 
241 {
243  scalarField& nu = tnu.ref();
244 
247 
248  nu = turbulenceModel().nu()().boundaryField()[patch_.index()];
249 
250  return tnu;
251 }
252 
253 
255 {
256  /*
257  const polyMesh& mesh = patch_.patch().boundaryMesh().mesh();
258  const compressible::turbulenceModel& turbulenceModel =
259  mesh.lookupObject<compressible::turbulenceModel>("turbulenceModel");
260  tmp<scalarField> talphaEff = turbulenceModel.alphaEff(patch_.index());
261  */
262 
263  tmp<scalarField> talphaEff(new scalarField(patch_.size(), Zero));
264 
266  << "no abstract thermalDiffusion is implemented. Returning zero field";
267 
268 
269  return talphaEff;
270 }
271 
272 
274 {
275  tmp<scalarField> twallDist(new scalarField(patch_.size(), Zero));
276  scalarField& wallDist = twallDist.ref();
277 
279 
280  return twallDist;
281 }
282 
283 
286 {
289 
290  tmp<scalarField> tdiffCoeff =
291  adjointRAS().diffusionCoeffVar1(patch_.index());
292 
293  return tdiffCoeff;
294 }
295 
296 
299 {
302 
303  tmp<scalarField> tdiffCoeff =
304  adjointRAS().diffusionCoeffVar2(patch_.index());
305 
306  return tdiffCoeff;
307 }
308 
309 
311 {
312  const autoPtr<incompressible::RASModelVariables>& RASVariables =
314  tmp<scalarField> tboundField(new scalarField(patch_.size(), Zero));
315  scalarField& boundField = tboundField.ref();
316 
317  boundField = RASVariables().TMVar1().boundaryField()[patch_.index()];
318 
319  return tboundField;
320 }
321 
322 
324 {
325  const autoPtr<incompressible::RASModelVariables>& RASVariables =
327  tmp<scalarField> tboundField(new scalarField(patch_.size(), Zero));
328  scalarField& boundField = tboundField.ref();
329 
330  boundField = RASVariables().TMVar2().boundaryField()[patch_.index()];
331 
332  return tboundField;
333 }
334 
335 
337 {
338  return primalVars_.U().boundaryField()[patch_.index()];
339 }
340 
341 
343 {
344  return primalVars_.p().boundaryField()[patch_.index()];
345 }
346 
347 
348 const fvsPatchScalarField&
350 {
351  return primalVars_.phi().boundaryField()[patch_.index()];
352 }
353 
354 
355 const fvPatchScalarField&
357 {
358  return
359  primalVars_.RASModelVariables()().nutRef().boundaryField()
360  [
361  patch_.index()
362  ];
363 }
364 
365 
367 {
368  return adjointVars().UaInst().boundaryField()[patch_.index()];
369 }
370 
371 
373 {
374  return adjointVars().paInst().boundaryField()[patch_.index()];
375 }
376 
377 
378 const fvsPatchScalarField&
380 {
382 }
383 
384 
386 {
387  return primalVars_.solverName();
388 }
389 
390 
392 {
393  return adjointVars().solverName();
394 }
395 
396 
397 const incompressibleVars&
399 {
400  return primalVars_;
401 }
402 
403 
406 {
408 }
409 
410 
413 {
414  return objectiveManager_;
415 }
416 
417 
418 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
419 
420 } // End namespace Foam
421 
422 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::fvPatchField< vector >
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::objectiveManager
class for managing incompressible objective functions.
Definition: objectiveManager.H:54
Foam::variablesSet::solverName
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
adjointRASModel.H
Foam::boundaryAdjointContributionIncompressible::laminarDiffusivity
tmp< scalarField > laminarDiffusivity()
Definition: boundaryAdjointContributionIncompressible.C:240
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::incompressibleAdjointSolver::getAdjointVars
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
Definition: incompressibleAdjointSolver.C:129
Foam::boundaryAdjointContributionIncompressible::velocitySource
tmp< vectorField > velocitySource()
Definition: boundaryAdjointContributionIncompressible.C:89
Foam::boundaryAdjointContributionIncompressible::phiab
const fvsPatchScalarField & phiab() const
Definition: boundaryAdjointContributionIncompressible.C:379
Foam::incompressibleAdjointMeanFlowVars::phiaInst
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
Definition: incompressibleAdjointMeanFlowVars.C:249
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::incompressibleVars::phi
const surfaceScalarField & phi() const
Return const reference to volume flux.
Definition: incompressibleVars.C:357
Foam::boundaryAdjointContributionIncompressible::adjointTMVariable1Source
tmp< scalarField > adjointTMVariable1Source()
Definition: boundaryAdjointContributionIncompressible.C:197
Foam::incompressibleVars::p
const volScalarField & p() const
Return const reference to pressure.
Definition: incompressibleVars.C:305
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::boundaryAdjointContributionIncompressible::normalVelocitySource
tmp< vectorField > normalVelocitySource()
Definition: boundaryAdjointContributionIncompressible.C:168
Foam::boundaryAdjointContributionIncompressible::Ub
const fvPatchVectorField & Ub() const
Definition: boundaryAdjointContributionIncompressible.C:336
boundaryAdjointContributionIncompressible.H
Foam::incompressibleVars::RASModelVariables
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Definition: incompressibleVars.C:444
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::fvPatch::size
virtual label size() const
Return size.
Definition: fvPatch.H:179
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:144
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::boundaryAdjointContributionIncompressible::energySource
tmp< scalarField > energySource()
Definition: boundaryAdjointContributionIncompressible.C:182
Foam::objectiveIncompressible::boundarydJdv
const boundaryVectorField & boundarydJdv()
Objective partial deriv wrt velocity for all patches.
Definition: objectiveIncompressible.C:378
Foam::boundaryAdjointContributionIncompressible::pab
const fvPatchScalarField & pab() const
Definition: boundaryAdjointContributionIncompressible.C:372
Foam::boundaryAdjointContributionIncompressible::adjointSolverName
const word adjointSolverName() const
Definition: boundaryAdjointContributionIncompressible.C:391
Foam::boundaryAdjointContributionIncompressible::wallDistance
tmp< scalarField > wallDistance()
Definition: boundaryAdjointContributionIncompressible.C:273
Foam::Field< vector >
Foam::wallDist
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition: wallDist.H:75
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
Foam::incompressibleAdjointMeanFlowVars::paInst
const volScalarField & paInst() const
Return const reference to pressure.
Definition: incompressibleAdjointMeanFlowVars.C:225
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::incompressibleVars::U
const volVectorField & U() const
Return const reference to velocity.
Definition: incompressibleVars.C:331
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::boundaryAdjointContributionIncompressible::phib
const fvsPatchScalarField & phib() const
Definition: boundaryAdjointContributionIncompressible.C:349
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::boundaryAdjointContributionIncompressible::Uab
const fvPatchVectorField & Uab() const
Definition: boundaryAdjointContributionIncompressible.C:366
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::boundaryAdjointContributionIncompressible::TMVariable1Diffusion
tmp< scalarField > TMVariable1Diffusion()
Definition: boundaryAdjointContributionIncompressible.C:285
Foam::objectiveIncompressible::boundarydJdT
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
Definition: objectiveIncompressible.C:418
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::boundaryAdjointContributionIncompressible::adjointTMVariable2Source
tmp< scalarField > adjointTMVariable2Source()
Definition: boundaryAdjointContributionIncompressible.C:212
Foam::incompressibleVars::turbulence
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Definition: incompressibleVars.C:431
Foam::objectiveIncompressible::boundarydJdvn
const boundaryScalarField & boundarydJdvn()
Objective partial deriv wrt normal velocity for all patches.
Definition: objectiveIncompressible.C:388
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::boundaryAdjointContribution
Abstract base class for computing contributions of the objective functions to the adjoint boundary co...
Definition: boundaryAdjointContribution.H:58
Foam::boundaryAdjointContributionIncompressible::adjointSolver_
const incompressibleAdjointSolver & adjointSolver_
Definition: boundaryAdjointContributionIncompressible.H:95
Foam::boundaryAdjointContributionIncompressible::momentumDiffusion
tmp< scalarField > momentumDiffusion()
Definition: boundaryAdjointContributionIncompressible.C:226
Foam::boundaryAdjointContributionIncompressible::TMVariable1
tmp< scalarField > TMVariable1()
Definition: boundaryAdjointContributionIncompressible.C:310
Foam::boundaryAdjointContribution::patch_
const fvPatch & patch_
Definition: boundaryAdjointContribution.H:78
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::boundaryAdjointContributionIncompressible::pb
const fvPatchScalarField & pb() const
Definition: boundaryAdjointContributionIncompressible.C:342
Foam::incompressibleAdjointVars::adjointTurbulence
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
Definition: incompressibleAdjointVars.C:70
Foam::boundaryAdjointContributionIncompressible::pressureSource
tmp< scalarField > pressureSource()
Definition: boundaryAdjointContributionIncompressible.C:110
Foam::objectiveIncompressible::boundarydJdTMvar2
const boundaryScalarField & boundarydJdTMvar2()
Objective partial deriv wrt turbulence model var 2 for all patches.
Definition: objectiveIncompressible.C:438
Foam::boundaryAdjointContributionIncompressible::TMVariable2
tmp< scalarField > TMVariable2()
Definition: boundaryAdjointContributionIncompressible.C:323
Foam::boundaryAdjointContributionIncompressible::adjointVars
const incompressibleAdjointVars & adjointVars() const
Definition: boundaryAdjointContributionIncompressible.C:405
Foam::boundaryAdjointContributionIncompressible::primalVars
const incompressibleVars & primalVars() const
Definition: boundaryAdjointContributionIncompressible.C:398
Foam::boundaryAdjointContributionIncompressible::sumContributions
tmp< Field< returnType > > sumContributions(PtrList< sourceType > &sourceList, const fvPatchField< returnType > &(castType::*boundaryFunction)(const label))
Foam::objectiveIncompressible::boundarydJdp
const boundaryVectorField & boundarydJdp()
Objective partial deriv wrt pressure (times normal) for all patches.
Definition: objectiveIncompressible.C:408
Foam::boundaryAdjointContributionIncompressible::tangentVelocitySource
tmp< vectorField > tangentVelocitySource()
Definition: boundaryAdjointContributionIncompressible.C:139
Foam::boundaryAdjointContributionIncompressible::turbulentDiffusivity
const fvPatchScalarField & turbulentDiffusivity() const
Definition: boundaryAdjointContributionIncompressible.C:356
Foam::objectiveIncompressible::boundarydJdTMvar1
const boundaryScalarField & boundarydJdTMvar1()
Objective partial deriv wrt turbulence model var 1 for all patches.
Definition: objectiveIncompressible.C:428
Foam::objectiveIncompressible::boundarydJdvt
const boundaryVectorField & boundarydJdvt()
Objective partial deriv wrt tangent velocity for all patches.
Definition: objectiveIncompressible.C:398
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::boundaryAdjointContributionIncompressible::getObjectiveManager
objectiveManager & getObjectiveManager()
Definition: boundaryAdjointContributionIncompressible.C:412
Foam::boundaryAdjointContributionIncompressible::objectiveManager_
objectiveManager & objectiveManager_
Definition: boundaryAdjointContributionIncompressible.H:86
Foam::boundaryAdjointContributionIncompressible::primalVars_
const incompressibleVars & primalVars_
Definition: boundaryAdjointContributionIncompressible.H:88
Foam::boundaryAdjointContributionIncompressible::TMVariable2Diffusion
tmp< scalarField > TMVariable2Diffusion()
Definition: boundaryAdjointContributionIncompressible.C:298
Foam::boundaryAdjointContributionIncompressible::thermalDiffusion
tmp< scalarField > thermalDiffusion()
Definition: boundaryAdjointContributionIncompressible.C:254
Foam::incompressibleAdjointMeanFlowVars::UaInst
const volVectorField & UaInst() const
Return const reference to velocity.
Definition: incompressibleAdjointMeanFlowVars.C:237
Foam::boundaryAdjointContributionIncompressible::primalSolverName
const word primalSolverName() const
Definition: boundaryAdjointContributionIncompressible.C:385
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::objectiveManager::getObjectiveFunctions
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Definition: objectiveManager.C:296