adjointSensitivityIncompressible.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 
30 #include "runTimeSelectionTables.H"
34 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 namespace incompressible
42 {
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(adjointSensitivity, 0);
47 defineRunTimeSelectionTable(adjointSensitivity, dictionary);
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 adjointSensitivity::adjointSensitivity
52 (
53  const fvMesh& mesh,
54  const dictionary& dict,
55  incompressibleVars& primalVars,
56  incompressibleAdjointVars& adjointVars,
59 )
60 :
62  derivatives_(0),
63  primalVars_(primalVars),
64  adjointVars_(adjointVars),
65  objectiveManager_(objectiveManager),
66  fvOptionsAdjoint_(fvOptionsAdjoint)
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
71 
73 (
74  const fvMesh& mesh,
75  const dictionary& dict,
76  incompressibleVars& primalVars,
77  incompressibleAdjointVars& adjointVars,
80 )
81 {
82  const word modelType(dict.get<word>("type"));
83 
84  Info<< "adjointSensitivity type : " << modelType << endl;
85 
86  auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
87 
88  if (!cstrIter.found())
89  {
91  (
92  dict,
93  "adjointSensitivity",
94  modelType,
95  *dictionaryConstructorTablePtr_
96  ) << exit(FatalIOError);
97  }
98 
100  (
101  cstrIter()
102  (
103  mesh,
104  dict,
105  primalVars,
106  adjointVars,
109  )
110  );
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
115 
117 {
119  write(type());
120  return derivatives_;
121 }
122 
123 
125 {
126  derivatives_ = scalar(0);
127  if (fieldSensPtr_.valid())
128  {
129  fieldSensPtr_().primitiveFieldRef() = scalar(0);
130  }
131 }
132 
133 
134 void adjointSensitivity::write(const word& baseName)
135 {
136  sensitivity::write(baseName);
137 }
138 
139 
141 {
142  // Term depending on the adjoint turbulence model
144  (
146  );
147  tmp<volTensorField> tturbulenceTerm(adjointRAS->FISensitivityTerm());
148  volTensorField& turbulenceTerm = tturbulenceTerm.ref();
149 
150  // nu effective
151  tmp<volScalarField> tnuEff(adjointRAS->nuEff());
152  const volScalarField& nuEff = tnuEff();
153 
154  tmp<volTensorField> tflowTerm
155  (
156  new volTensorField
157  (
158  IOobject
159  (
160  "flowTerm",
161  mesh_.time().timeName(),
162  mesh_,
165  ),
166  mesh_,
168  )
169  );
170  volTensorField& flowTerm = tflowTerm.ref();
171 
172  const volScalarField& p = primalVars_.p();
173  const volVectorField& U = primalVars_.U();
174  const volScalarField& pa = adjointVars_.pa();
175  const volVectorField& Ua = adjointVars_.Ua();
176  volTensorField gradU(fvc::grad(U));
177  volTensorField gradUa(fvc::grad(Ua));
178 
179  // Explicitly correct the boundary gradient to get rid of
180  // the tangential component
181  forAll(mesh_.boundary(), patchI)
182  {
183  const fvPatch& patch = mesh_.boundary()[patchI];
184  if (isA<wallFvPatch>(patch))
185  {
186  tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
187  const vectorField& nf = tnf();
188  gradU.boundaryFieldRef()[patchI] =
189  nf*U.boundaryField()[patchI].snGrad();
190  //gradUa.boundaryField()[patchI] =
191  // nf*Ua.boundaryField()[patchI].snGrad();
192  }
193  }
194 
195  volTensorField stress(nuEff*(gradU + T(gradU)));
196  autoPtr<volVectorField> stressXPtr
197  (
198  createZeroFieldPtr<vector>(mesh_, "stressX", stress.dimensions())
199  );
200  autoPtr<volVectorField> stressYPtr
201  (
202  createZeroFieldPtr<vector>(mesh_, "stressY", stress.dimensions())
203  );
204  autoPtr<volVectorField> stressZPtr
205  (
206  createZeroFieldPtr<vector>(mesh_, "stressZ", stress.dimensions())
207  );
208 
209  stressXPtr().replace(0, stress.component(0));
210  stressXPtr().replace(1, stress.component(1));
211  stressXPtr().replace(2, stress.component(2));
212 
213  stressYPtr().replace(0, stress.component(3));
214  stressYPtr().replace(1, stress.component(4));
215  stressYPtr().replace(2, stress.component(5));
216 
217  stressZPtr().replace(0, stress.component(6));
218  stressZPtr().replace(1, stress.component(7));
219  stressZPtr().replace(2, stress.component(8));
220 
221  volTensorField gradStressX(fvc::grad(stressXPtr()));
222  volTensorField gradStressY(fvc::grad(stressYPtr()));
223  volTensorField gradStressZ(fvc::grad(stressZPtr()));
224 
225  // Contribution from objective functions and constraints
226  volTensorField objectiveContributions
227  (
228  IOobject
229  (
230  "objectiveContributions",
231  mesh_.time().timeName(),
232  mesh_,
235  ),
236  mesh_,
238  );
240  forAll(functions, funcI)
241  {
242  objectiveContributions +=
243  functions[funcI].weight()
244  *functions[funcI].gradDxDbMultiplier();
245  }
246 
247  // Note:
248  // term4 (Ua & grad(stress)) is numerically tricky. Its div leads to third
249  // order spatial derivs in E-SI based computations Applying the product
250  // derivative rule (putting Ua inside the grad) gives better results in
251  // NACA0012, SA, WF. However, the original formulation should be kept at
252  // the boundary in order to respect the Ua boundary conditions (necessary
253  // for E-SI to give the same sens as FI). A mixed approach is hence
254  // followed
255  volTensorField term4
256  (
257  - nuEff*(gradUa & (gradU + T(gradU)))
258  + fvc::grad(nuEff * Ua & (gradU + T(gradU)))
259  );
260 
261  forAll(mesh_.boundary(), pI)
262  {
263  if (!isA<coupledFvPatch>(mesh_.boundary()[pI]))
264  {
265  term4.boundaryFieldRef()[pI] =
266  Ua.component(0)().boundaryField()[pI]
267  *gradStressX.boundaryField()[pI]
268  + Ua.component(1)().boundaryField()[pI]
269  *gradStressY.boundaryField()[pI]
270  + Ua.component(2)().boundaryField()[pI]
271  *gradStressZ.boundaryField()[pI];
272  }
273  }
274 
275  const autoPtr<ATCModel>& ATCModel =
277  (
279  ).getATCModel();
280 
281  // Compute dxdb multiplier
282  flowTerm =
283  // Term 1, ATC
285  // Term 2
286  - fvc::grad(p) * Ua
287  // Term 3
288  - nuEff*(gradU & (gradUa + T(gradUa)))
289  // Term 4
290  + term4
291  // Term 5
292  + (pa * gradU)
293  // Term 6, from the adjoint turbulence model
294  + turbulenceTerm.T()
295  // Term 7, term from objective functions
296  + objectiveContributions;
297 
298  flowTerm.correctBoundaryConditions();
299 
300  return (tflowTerm);
301 }
302 
303 
305 {
307  volTensorField& gradDxDbMult = tgradDxDbMult.ref();
308 
309  tmp<volVectorField> tadjointMeshMovementSource
310  (
311  new volVectorField
312  (
313  IOobject
314  (
315  "adjointMeshMovementSource",
316  mesh_.time().timeName(),
317  mesh_,
320  ),
321  mesh_,
322  dimensionedVector(gradDxDbMult.dimensions()/dimLength, Zero)
323  )
324  );
325 
326  volVectorField& source = tadjointMeshMovementSource.ref();
327 
328  source -= fvc::div(gradDxDbMult.T());
329 
330  // Terms from fvOptions
332  {
333  source += fvOptionsAdjoint_[oI].dxdbMult(adjointVars_);
334  }
335 
336  return (tadjointMeshMovementSource);
337 }
338 
339 
340 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 
342 } // End namespace incompressible
343 } // End namespace Foam
344 
345 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::incompressible::adjointSensitivity::derivatives_
scalarField derivatives_
Definition: adjointSensitivityIncompressible.H:84
Foam::objectiveManager
class for managing incompressible objective functions.
Definition: objectiveManager.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::incompressible::adjointSensitivity::calculateSensitivities
virtual const scalarField & calculateSensitivities()
Calculates and returns sensitivity fields.
Definition: adjointSensitivityIncompressible.C:116
Foam::incompressibleAdjointSolver
Base class for incompressibleAdjoint solvers.
Definition: incompressibleAdjointSolver.H:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::incompressible::adjointSensitivity::fvOptionsAdjoint_
fv::optionAdjointList & fvOptionsAdjoint_
Definition: adjointSensitivityIncompressible.H:88
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
boundaryAdjointContribution.H
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::incompressible::adjointSensitivity::assembleSensitivities
virtual void assembleSensitivities()=0
Assemble sensitivities.
wallFvPatch.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::incompressibleVars::p
const volScalarField & p() const
Return const reference to pressure.
Definition: incompressibleVars.C:305
Foam::incompressible::defineRunTimeSelectionTable
defineRunTimeSelectionTable(adjointSensitivity, dictionary)
Foam::ATCModel
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition: ATCModel.H:60
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::fv::optionAdjointList
Definition: fvOptionAdjointList.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
incompressibleAdjointSolver.H
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:380
Foam::sensitivity::mesh_
const fvMesh & mesh_
Definition: sensitivity.H:69
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::ATCModel::getFISensitivityTerm
virtual tmp< volTensorField > getFISensitivityTerm() const =0
Get the FI sensitivity derivatives term coming from the ATC.
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::incompressible::adjointSensitivity::New
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, incompressibleVars &primalVars, incompressibleAdjointVars &adjointVars, objectiveManager &objectiveManager, fv::optionAdjointList &fvOptionsAdjoint)
Return a reference to the selected turbulence model.
Definition: adjointSensitivityIncompressible.C:73
Foam::incompressibleAdjointMeanFlowVars::pa
const volScalarField & pa() const
Return const reference to pressure.
Definition: incompressibleAdjointMeanFlowVars.C:147
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:65
Foam::incompressible::adjointSensitivity::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: adjointSensitivityIncompressible.C:124
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::incompressibleAdjointMeanFlowVars::Ua
const volVectorField & Ua() const
Return const reference to velocity.
Definition: incompressibleAdjointMeanFlowVars.C:173
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::sensitivity::fieldSensPtr_
autoPtr< volScalarField > fieldSensPtr_
Definition: sensitivity.H:74
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::incompressible::adjointSensitivity::adjointVars_
incompressibleAdjointVars & adjointVars_
Definition: adjointSensitivityIncompressible.H:86
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressible::adjointSensitivity::primalVars_
incompressibleVars & primalVars_
Definition: adjointSensitivityIncompressible.H:85
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::T
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Definition: GeometricField.C:1015
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
adjointSensitivityIncompressible.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:718
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::incompressibleAdjointVars::adjointTurbulence
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
Definition: incompressibleAdjointVars.C:70
Foam::incompressible::adjointSensitivity::adjointMeshMovementSource
tmp< volVectorField > adjointMeshMovementSource()
Compute source term for adjoint mesh movement equation.
Definition: adjointSensitivityIncompressible.C:304
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::sensitivity::write
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
Definition: sensitivity.C:77
Foam::objectiveManager::adjointSolverName
const word & adjointSolverName() const
Return name of the adjointSolver.
Definition: objectiveManager.C:255
Foam::GeometricField< tensor, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::sensitivity
Abstract base class for adjoint sensitivities.
Definition: sensitivity.H:63
fvOptionsAdjoint
fv::IOoptionListAdjoint fvOptionsAdjoint(mesh)
Foam::dimensionedTensor
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedTensor.H:51
Foam::incompressible::adjointSensitivity::computeGradDxDbMultiplier
tmp< volTensorField > computeGradDxDbMultiplier()
Definition: adjointSensitivityIncompressible.C:140
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::incompressible::adjointSensitivity::write
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
Definition: adjointSensitivityIncompressible.C:134
Foam::objectiveManager::getObjectiveFunctions
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Definition: objectiveManager.C:243
Foam::incompressible::adjointSensitivity::objectiveManager_
objectiveManager & objectiveManager_
Definition: adjointSensitivityIncompressible.H:87