Loading...
Searching...
No Matches
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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 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#include "adjointRASModel.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
52/*---------------------------------------------------------------------------*\
53 Class objectiveIncompressible Declaration
54\*---------------------------------------------------------------------------*/
55
56class objectiveIncompressible
57:
58 public objective
60protected:
61
62 // Protected data
63
66 // Contribution to field adjoint equations
67 // v,p,T and turbulence model variables
68 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
72
73 //- First turbulence model variable
75
76 //- Second turbulence model variable
78
79 // Contribution to adjoint boundary conditions
80 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
82
83 //- Adjoint outlet pressure
86 //- Adjoint outlet velocity
88
89 //- Adjoint (intlet,wall) velocity
91
92 //- Adjoint outlet temperature
94
95 //- Adjoint outlet turbulence model var 1
97
98 //- Adjoint outlet turbulence model var 2
101 //- Jacobian wrt to nut
103
104 //- Term multiplying gradU variations
106
107
108private:
109
110 // Private Member Functions
111
112 //- No copy construct
113 objectiveIncompressible(const objectiveIncompressible&) = delete;
114
115 //- No copy assignment
116 void operator=(const objectiveIncompressible&) = delete;
117
118
119public:
121 //- Runtime type information
122 TypeName("incompressible");
123
124
125 // Declare run-time constructor selection table
126
128 (
129 autoPtr,
130 objectiveIncompressible,
132 (
133 const fvMesh& mesh,
134 const dictionary& dict,
135 const word& adjointSolverName,
136 const word& primalSolverName
137 ),
138 (mesh, dict, adjointSolverName, primalSolverName)
139 );
140
141
142 // Constructors
144 //- Construct from components
145 objectiveIncompressible
146 (
147 const fvMesh& mesh,
149 const word& adjointSolverName,
150 const word& primalSolverName
151 );
152
153
154 // Selectors
155
156 //- Return a reference to the selected turbulence model
158 (
159 const fvMesh& mesh,
160 const dictionary& dict,
161 const word& adjointSolverName,
162 const word& primalSolverName
163 );
164
165
166 //- Destructor
167 virtual ~objectiveIncompressible() = default;
168
169
170 // Member Functions
171
172 //- Return the objective function value
173 virtual scalar J() = 0;
174
175 //- Normalize all fields allocated by the objective
176 virtual void doNormalization();
177
178 //- Contribution to field adjoint momentum eqs
179 inline const volVectorField& dJdv();
180
181 //- Contribution to field adjoint continuity eq
182 inline const volScalarField& dJdp();
183
184 //- Contribution to field adjoint energy eq
185 inline const volScalarField& dJdT();
186
187 //- Contribution to field adjoint turbulence model variable 1
188 inline const volScalarField& dJdTMvar1();
189
190 //- Contribution to field adjoint turbulence model variable 2
191 inline const volScalarField& dJdTMvar2();
192
193 //- Objective partial deriv wrt velocity for a specific patch
194 const fvPatchVectorField& boundarydJdv(const label);
195
196 //- Objective partial deriv wrt normal velocity for a specific patch
197 inline const fvPatchScalarField& boundarydJdvn(const label);
198
199 //- Objective partial deriv wrt tangent velocity for a specific patch
200 inline const fvPatchVectorField& boundarydJdvt(const label);
201
202 //- Objective partial deriv wrt pressure (times normal) for a specific
203 //- patch
204 inline const fvPatchVectorField& boundarydJdp(const label);
205
206 //- Objective partial deriv wrt temperature for a specific patch
207 inline const fvPatchScalarField& boundarydJdT(const label);
208
209 //- Objective partial deriv wrt turbulence model var 1 for a specific
210 //- patch
211 inline const fvPatchScalarField& boundarydJdTMvar1(const label);
212
213 //- Objective partial deriv wrt turbulence model var 2 for a specific
214 //- patch
215 inline const fvPatchScalarField& boundarydJdTMvar2(const label);
216
217 //- Objective partial deriv wrt nut for a specific patch
218 inline const fvPatchScalarField& boundarydJdnut(const label);
219
220 //- Objective partial deriv wrt stress tensor
221 inline const fvPatchTensorField& boundarydJdGradU(const label);
222
223 //- Objective partial deriv wrt velocity for all patches
224 inline const boundaryVectorField& boundarydJdv();
225
226 //- Objective partial deriv wrt normal velocity for all patches
227 inline const boundaryScalarField& boundarydJdvn();
228
229 //- Objective partial deriv wrt tangent velocity for all patches
230 inline const boundaryVectorField& boundarydJdvt();
231
232 //- Objective partial deriv wrt pressure (times normal) for all patches
233 inline const boundaryVectorField& boundarydJdp();
234
235 //- Objective partial deriv wrt temperature for all patches
236 inline const boundaryScalarField& boundarydJdT();
237
238 //- Objective partial deriv wrt turbulence model var 1 for all patches
240
241 //- Objective partial deriv wrt turbulence model var 2 for all patches
243
244 //- Objective partial deriv wrt nut for all patches
245 inline const boundaryScalarField& boundarydJdnut();
246
247 //- Objective partial deriv wrt gradU
249
250 //- Update objective function derivatives
251 virtual void update();
252
253 //- Update objective function derivatives
254 virtual void nullify();
255
256 // Funtions used in objectives including turbulent quantities
257
258 //- Allocate fields related to the differentiation of turbulence
259 //- models, if necessary
261
262 //- Check if cellZones provided include at least one cell
263 void checkCellZonesSize(const labelList& zoneIDs) const;
264
265 //- Compute dJdTMVar{1,2}
266 void update_dJdTMvar
267 (
268 autoPtr<volScalarField>& dJdTMvarPtr,
271 const,
272 const volScalarField& JacobianMultiplier,
273 const labelList& zones
274 );
275
276
277 //- Update vol and boundary fields and derivatives
278 // Do nothing in the base. The relevant ones should be overwritten
279 // in the child objective functions
280 virtual void update_dJdv()
281 {}
282
283 virtual void update_dJdp()
284 {}
285
286 virtual void update_dJdT()
287 {}
288
289 virtual void update_dJdTMvar1()
290 {}
291
292 virtual void update_dJdTMvar2()
293 {}
294
295 virtual void update_dJdb()
296 {}
297
298 virtual void update_dJdbField()
299 {}
300
301 virtual void update_divDxDbMultiplier()
302 {}
303
304 virtual void update_gradDxDbMultiplier()
305 {}
306
307 virtual void update_boundarydJdv()
308 {}
309
310 virtual void update_boundarydJdvn()
311 {}
312
313 virtual void update_boundarydJdvt()
314 {}
315
316 virtual void update_boundarydJdp()
317 {}
318
319 virtual void update_boundarydJdT()
320 {}
321
322 virtual void update_boundarydJdTMvar1()
323 {}
324
325 virtual void update_boundarydJdTMvar2()
326 {}
327
328 virtual void update_boundarydJdnut()
329 {}
330
331 virtual void update_boundarydJdGradU()
332 {}
333
334 //- Vector sources can be given only to the adjoint momentum equations.
335 //- Implemented in base objectiveIncompressible
336 virtual void addSource(fvVectorMatrix& matrix);
337
338 //- Scalar sources are more ambigious (adjoint pressure, turbulence
339 //- model, energy, etc equations), so the equivalent functions should
340 //- be overridden on an objective-basis
341 virtual void addSource(fvScalarMatrix& matrix)
342 {}
343
344
345 //- Some objectives need to store some auxiliary values.
346 //- If averaging is enabled, update these mean values here.
347 // By convention, the mean values (eg mean drag) refer to these flow
348 // values computed using the mean fields, rather than averaging the
349 // values themselves
350 virtual void update_meanValues()
351 {}
352
353 //- Write objective function history
354 virtual bool write(const bool valid = true) const;
355
356 //- Inline functions for checking whether pointers are set or not
357 inline bool hasdJdv() const noexcept;
358 inline bool hasdJdp() const noexcept;
359 inline bool hasdJdT() const noexcept;
360 inline bool hasdJdTMVar1() const noexcept;
361 inline bool hasdJdTMVar2() const noexcept;
362 inline bool hasBoundarydJdv() const noexcept;
363 inline bool hasBoundarydJdvn() const noexcept;
364 inline bool hasBoundarydJdvt() const noexcept;
365 inline bool hasBoundarydJdp() const noexcept;
366 inline bool hasBoundarydJdT() const noexcept;
367 inline bool hasBoundarydJdTMVar1() const noexcept;
368 inline bool hasBoundarydJdTMVar2() const noexcept;
369 inline bool hasBoundarydJdnut() const noexcept;
370 inline bool hasBoundarydJdGradU() const noexcept;
371};
372
374// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375
376} // End namespace Foam
377
378// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
380#include "objectiveIncompressibleI.H"
381
382// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
383
384#endif
386// ************************************************************************* //
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for incompressible turbulence models.
Base class for solution control classes.
const volScalarField & dJdTMvar1()
Contribution to field adjoint turbulence model variable 1.
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
const volScalarField & dJdT()
Contribution to field adjoint energy eq.
virtual void addSource(fvScalarMatrix &matrix)
Scalar sources are more ambigious (adjoint pressure, turbulence model, energy, etc equations),...
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
static autoPtr< objectiveIncompressible > New(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
bool hasdJdv() const noexcept
Inline functions for checking whether pointers are set or not.
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
TypeName("incompressible")
Runtime type information.
autoPtr< boundaryScalarField > bdJdTMvar2Ptr_
Adjoint outlet turbulence model var 2.
autoPtr< volScalarField > dJdpPtr_
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
virtual void addSource(fvVectorMatrix &matrix)
Vector sources can be given only to the adjoint momentum equations. Implemented in base objectiveInco...
const volScalarField & dJdTMvar2()
Contribution to field adjoint turbulence model variable 2.
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
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.
const volVectorField & dJdv()
Contribution to field adjoint momentum eqs.
declareRunTimeSelectionTable(autoPtr, objectiveIncompressible, dictionary,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
virtual void update_meanValues()
Some objectives need to store some auxiliary values. If averaging is enabled, update these mean value...
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
const volScalarField & dJdp()
Contribution to field adjoint continuity eq.
autoPtr< boundaryScalarField > bdJdTMvar1Ptr_
Adjoint outlet turbulence model var 1.
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
const boundaryVectorField & boundarydJdp()
Objective partial deriv wrt pressure (times normal) for all patches.
const boundaryVectorField & boundarydJdvt()
Objective partial deriv wrt tangent velocity for all patches.
virtual scalar J()=0
Return the objective function value.
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
autoPtr< boundaryScalarField > bdJdTPtr_
Adjoint outlet temperature.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
autoPtr< volScalarField > dJdTPtr_
const boundaryTensorField & boundarydJdGradU()
Objective partial deriv wrt gradU.
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
virtual void doNormalization()
Normalize all fields allocated by the objective.
const boundaryScalarField & boundarydJdnut()
Objective partial deriv wrt nut for all patches.
virtual void update()
Update objective function derivatives.
const boundaryScalarField & boundarydJdTMvar1()
Objective partial deriv wrt turbulence model var 1 for all patches.
const boundaryScalarField & boundarydJdTMvar2()
Objective partial deriv wrt turbulence model var 2 for all patches.
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
const incompressibleVars & vars_
autoPtr< boundaryVectorField > bdJdvPtr_
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
const boundaryScalarField & boundarydJdvn()
Objective partial deriv wrt normal velocity for all patches.
autoPtr< volVectorField > dJdvPtr_
const dictionary & dict() const
Return objective dictionary.
Definition objective.C:91
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
const labelIOList & zoneIDs
Definition correctPhi.H:59
Namespace for OpenFOAM.
volScalarField::Boundary boundaryScalarField
volFields
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fvPatchField< tensor > fvPatchTensorField
volTensorField::Boundary boundaryTensorField
fvMatrix< vector > fvVectorMatrix
const direction noexcept
Definition scalarImpl.H:265
volVectorField::Boundary boundaryVectorField
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
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:68