objectiveIncompressible.H
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-2020 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
29Class
30 Foam::objectiveIncompressible
31
32Description
33 Abstract base class for objective functions in incompressible flows
34
35SourceFiles
36 objectiveIncompressible.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef objectiveIncompressible_H
41#define objectiveIncompressible_H
42
43#include "objective.H"
44#include "incompressibleVars.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class objectiveIncompressible Declaration
53\*---------------------------------------------------------------------------*/
56:
57 public objective
58{
59protected:
60
61 // Protected data
64
65 // Contribution to field adjoint equations
66 // v,p,T and turbulence model variables
67 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
71
72 //- First turbulence model variable
74
75 //- Second turbulence model variable
77
78 // Contribution to adjoint boundary conditions
79 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
81
82 //- Adjoint outlet pressure
84
85 //- Adjoint outlet velocity
87
88 //- Adjoint (intlet,wall) velocity
90
91 //- Adjoint outlet temperature
93
94 //- Adjoint outlet turbulence model var 1
96
97 //- Adjoint outlet turbulence model var 2
99
100 //- Jacobian wrt to nut
102
103 //- Term multiplying gradU variations
105
106
107private:
108
109 // Private Member Functions
110
111 //- No copy construct
113
114 //- No copy assignment
115 void operator=(const objectiveIncompressible&) = delete;
116
117
118public:
119
120 //- Runtime type information
121 TypeName("incompressible");
122
123
124 // Declare run-time constructor selection table
127 (
128 autoPtr,
131 (
132 const fvMesh& mesh,
133 const dictionary& dict,
134 const word& adjointSolverName,
135 const word& primalSolverName
136 ),
137 (mesh, dict, adjointSolverName, primalSolverName)
138 );
139
140
141 // Constructors
142
143 //- Construct from components
145 (
146 const fvMesh& mesh,
147 const dictionary& dict,
148 const word& adjointSolverName,
149 const word& primalSolverName
150 );
151
152
153 // Selectors
154
155 //- Return a reference to the selected turbulence model
157 (
158 const fvMesh& mesh,
159 const dictionary& dict,
160 const word& adjointSolverName,
161 const word& primalSolverName
162 );
163
164
165 //- Destructor
166 virtual ~objectiveIncompressible() = default;
167
168
169 // Member Functions
170
171 //- Return the objective function value
172 virtual scalar J() = 0;
173
174 //- Normalize all fields allocated by the objective
175 virtual void doNormalization();
176
177 //- Contribution to field adjoint momentum eqs
178 const volVectorField& dJdv();
179
180 //- Contribution to field adjoint continuity eq
181 const volScalarField& dJdp();
182
183 //- Contribution to field adjoint energy eq
184 const volScalarField& dJdT();
185
186 //- Contribution to field adjoint turbulence model variable 1
187 const volScalarField& dJdTMvar1();
188
189 //- Contribution to field adjoint turbulence model variable 2
190 const volScalarField& dJdTMvar2();
191
192 //- Objective partial deriv wrt velocity for a specific patch
193 const fvPatchVectorField& boundarydJdv(const label);
194
195 //- Objective partial deriv wrt normal velocity for a specific patch
196 const fvPatchScalarField& boundarydJdvn(const label);
197
198 //- Objective partial deriv wrt tangent velocity for a specific patch
199 const fvPatchVectorField& boundarydJdvt(const label);
200
201 //- Objective partial deriv wrt pressure (times normal) for a specific
202 //- patch
203 const fvPatchVectorField& boundarydJdp(const label);
204
205 //- Objective partial deriv wrt temperature for a specific patch
206 const fvPatchScalarField& boundarydJdT(const label);
207
208 //- Objective partial deriv wrt turbulence model var 1 for a specific
209 //- patch
210 const fvPatchScalarField& boundarydJdTMvar1(const label);
211
212 //- Objective partial deriv wrt turbulence model var 2 for a specific
213 //- patch
214 const fvPatchScalarField& boundarydJdTMvar2(const label);
215
216 //- Objective partial deriv wrt nut for a specific patch
217 const fvPatchScalarField& boundarydJdnut(const label);
218
219 //- Objective partial deriv wrt stress tensor
220 const fvPatchTensorField& boundarydJdGradU(const label);
221
222 //- Objective partial deriv wrt velocity for all patches
224
225 //- Objective partial deriv wrt normal velocity for all patches
227
228 //- Objective partial deriv wrt tangent velocity for all patches
230
231 //- Objective partial deriv wrt pressure (times normal) for all patches
233
234 //- Objective partial deriv wrt temperature for all patches
236
237 //- Objective partial deriv wrt turbulence model var 1 for all patches
239
240 //- Objective partial deriv wrt turbulence model var 2 for all patches
242
243 //- Objective partial deriv wrt nut for all patches
245
246 //- Objective partial deriv wrt gradU
248
249 //- Update objective function derivatives
250 virtual void update();
251
252 //- Update objective function derivatives
253 virtual void nullify();
254
255 //- Update vol and boundary fields and derivatives
256 // Do nothing in the base. The relevant ones should be overwritten
257 // in the child objective functions
258 virtual void update_dJdv()
259 {}
261 virtual void update_dJdp()
262 {}
264 virtual void update_dJdT()
265 {}
267 virtual void update_dJdTMvar1()
268 {}
270 virtual void update_dJdTMvar2()
271 {}
273 virtual void update_dJdb()
274 {}
276 virtual void update_divDxDbMultiplier()
277 {}
279 virtual void update_gradDxDbMultiplier()
280 {}
282 virtual void update_boundarydJdv()
283 {}
285 virtual void update_boundarydJdvn()
286 {}
288 virtual void update_boundarydJdvt()
289 {}
291 virtual void update_boundarydJdp()
292 {}
294 virtual void update_boundarydJdT()
295 {}
297 virtual void update_boundarydJdTMvar1()
298 {}
300 virtual void update_boundarydJdTMvar2()
301 {}
303 virtual void update_boundarydJdnut()
304 {}
306 virtual void update_boundarydJdGradU()
307 {}
309 virtual void update_boundarydJdb()
310 {}
312 virtual void update_dSdbMultiplier()
313 {}
315 virtual void update_dndbMultiplier()
316 {}
318 virtual void update_dxdbMultiplier()
319 {}
321 virtual void update_dxdbDirectMultiplier()
322 {}
323
324 //- Some objectives need to store some auxiliary values.
325 //- If averaging is enabled, update these mean values here.
326 // By convention, the mean values (eg mean drag) refer to these flow
327 // values computed using the mean fields, rather than averaging the
328 // values themselves
329 virtual void update_meanValues()
330 {}
331
332 //- Write objective function history
333 virtual bool write(const bool valid = true) const;
334
335 //- Inline functions for checking whether pointers are set or not
336 inline bool hasdJdv() const;
337 inline bool hasdJdp() const;
338 inline bool hasdJdT() const;
339 inline bool hasdJdTMVar1() const;
340 inline bool hasdJdTMVar2() const;
341 inline bool hasBoundarydJdv() const;
342 inline bool hasBoundarydJdvn() const;
343 inline bool hasBoundarydJdvt() const;
344 inline bool hasBoundarydJdp() const;
345 inline bool hasBoundarydJdT() const;
346 inline bool hasBoundarydJdTMVar1() const;
347 inline bool hasBoundarydJdTMVar2() const;
348 inline bool hasBoundarydJdnut() const;
349 inline bool hasBoundarydJdGradU() const;
350};
351
352
353// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354
355} // End namespace Foam
356
357// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
358
360
361// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362
363#endif
364
365// ************************************************************************* //
Generic GeometricBoundaryField class.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for solution control classes.
Abstract base class for objective functions in incompressible flows.
const boundaryVectorField & boundarydJdp()
Objective partial deriv wrt pressure (times normal) for all patches.
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
const volScalarField & dJdTMvar1()
Contribution to field adjoint turbulence model variable 1.
const boundaryScalarField & boundarydJdTMvar1()
Objective partial deriv wrt turbulence model var 1 for all patches.
static autoPtr< objectiveIncompressible > New(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
TypeName("incompressible")
Runtime type information.
autoPtr< boundaryScalarField > bdJdTMvar2Ptr_
Adjoint outlet turbulence model var 2.
autoPtr< volScalarField > dJdpPtr_
const boundaryTensorField & boundarydJdGradU()
Objective partial deriv wrt gradU.
const volScalarField & dJdT()
Contribution to field adjoint energy eq.
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
virtual void update_boundarydJdb()
Update objective function derivative term.
const boundaryScalarField & boundarydJdnut()
Objective partial deriv wrt nut for all patches.
virtual void update_dJdv()
Update vol and boundary fields and derivatives.
virtual void nullify()
Update objective function derivatives.
const boundaryVectorField & boundarydJdv()
Objective partial deriv wrt velocity for all patches.
virtual ~objectiveIncompressible()=default
Destructor.
declareRunTimeSelectionTable(autoPtr, objectiveIncompressible, dictionary,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
autoPtr< boundaryScalarField > bdJdTMvar1Ptr_
Adjoint outlet turbulence model var 1.
const boundaryScalarField & boundarydJdvn()
Objective partial deriv wrt normal velocity for all patches.
virtual scalar J()=0
Return the objective function value.
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
autoPtr< boundaryScalarField > bdJdTPtr_
Adjoint outlet temperature.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
autoPtr< volScalarField > dJdTPtr_
const boundaryVectorField & boundarydJdvt()
Objective partial deriv wrt tangent velocity for all patches.
const volScalarField & dJdp()
Contribution to field adjoint continuity eq.
const volScalarField & dJdTMvar2()
Contribution to field adjoint turbulence model variable 2.
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
virtual void doNormalization()
Normalize all fields allocated by the objective.
virtual void update()
Update objective function derivatives.
bool hasdJdv() const
Inline functions for checking whether pointers are set or not.
const volVectorField & dJdv()
Contribution to field adjoint momentum eqs.
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
const incompressibleVars & vars_
const boundaryScalarField & boundarydJdTMvar2()
Objective partial deriv wrt turbulence model var 2 for all patches.
autoPtr< boundaryVectorField > bdJdvPtr_
autoPtr< volVectorField > dJdvPtr_
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:94
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Namespace for OpenFOAM.
runTime write()
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73