objectiveIncompressible.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-2021 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\*---------------------------------------------------------------------------*/
29
32#include "createZeroField.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
45(
49);
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
54(
55 const fvMesh& mesh,
56 const dictionary& dict,
57 const word& adjointSolverName,
58 const word& primalSolverName
59)
60:
61 objective(mesh, dict, adjointSolverName, primalSolverName),
62
63 vars_
64 (
65 mesh.lookupObject<incompressiblePrimalSolver>(primalSolverName).
66 getIncoVars()
67 ),
68
69 // Initialize pointers to nullptr.
70 // Not all of them are required for each objective function.
71 // Each child should allocate whatever is needed.
72
73 // Field adjoint Eqs
74 dJdvPtr_(nullptr),
75 dJdpPtr_(nullptr),
76 dJdTPtr_(nullptr),
77 dJdTMvar1Ptr_(nullptr),
78 dJdTMvar2Ptr_(nullptr),
79
80 // Adjoint boundary conditions
81 bdJdvPtr_(nullptr),
82 bdJdvnPtr_(nullptr),
83 bdJdvtPtr_(nullptr),
84 bdJdpPtr_(nullptr),
85 bdJdTPtr_(nullptr),
86 bdJdTMvar1Ptr_(nullptr),
87 bdJdTMvar2Ptr_(nullptr),
88 bdJdnutPtr_(nullptr),
89 bdJdGradUPtr_(nullptr)
90{
91 weight_ = dict.get<scalar>("weight");
93}
94
95
96// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
97
99(
100 const fvMesh& mesh,
101 const dictionary& dict,
102 const word& adjointSolverName,
103 const word& primalSolverName
104)
105{
106 const word modelType(dict.get<word>("type"));
107
108 Info<< "Creating objective function : " << dict.dictName()
109 << " of type " << modelType << endl;
110
111 auto* ctorPtr = dictionaryConstructorTable(modelType);
112
113 if (!ctorPtr)
114 {
116 (
117 dict,
118 "objectiveIncompressible",
119 modelType,
120 *dictionaryConstructorTablePtr_
121 ) << exit(FatalIOError);
122 }
123
125 (
126 ctorPtr(mesh, dict, adjointSolverName, primalSolverName)
127 );
128}
129
130
131// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132
134{
135 if (normalize_ && normFactor_)
136 {
137 const scalar oneOverNorm(1./normFactor_());
138
139 if (hasdJdv())
140 {
141 dJdvPtr_().primitiveFieldRef() *= oneOverNorm;
142 }
143 if (hasdJdp())
144 {
145 dJdpPtr_().primitiveFieldRef() *= oneOverNorm;
146 }
147 if (hasdJdT())
148 {
149 dJdTPtr_().primitiveFieldRef() *= oneOverNorm;
150 }
151 if (hasdJdTMVar1())
152 {
153 dJdTMvar1Ptr_().primitiveFieldRef() *= oneOverNorm;
154 }
155 if (hasdJdTMVar2())
156 {
157 dJdTMvar2Ptr_().primitiveFieldRef() *= oneOverNorm;
158 }
159 if (hasBoundarydJdv())
160 {
161 bdJdvPtr_() *= oneOverNorm;
162 }
163 if (hasBoundarydJdvn())
164 {
165 bdJdvnPtr_() *= oneOverNorm;
166 }
167 if (hasBoundarydJdvt())
168 {
169 bdJdvtPtr_() *= oneOverNorm;
170 }
171 if (hasBoundarydJdp())
172 {
173 bdJdpPtr_() *= oneOverNorm;
174 }
175 if (hasBoundarydJdT())
176 {
177 bdJdTPtr_() *= oneOverNorm;
178 }
180 {
181 bdJdTMvar1Ptr_() *= oneOverNorm;
182 }
184 {
185 bdJdTMvar2Ptr_() *= oneOverNorm;
186 }
187 if (hasBoundarydJdnut())
188 {
189 bdJdnutPtr_() *= oneOverNorm;
190 }
192 {
193 bdJdGradUPtr_() *= oneOverNorm;
194 }
195
196 // Normalize geometric fields
198 }
199}
200
201
203{
204 if (!dJdvPtr_)
205 {
206 // If pointer is not set, set it to a zero field
207 dJdvPtr_.reset
208 (
209 createZeroFieldPtr<vector>
210 (
211 mesh_,
212 ("dJdv_"+type()),
213 dimensionSet(0, 3, -2, 0, 0, 0, 0)
214 )
215 );
216 }
217 return *dJdvPtr_;
218}
219
220
222{
223 if (!dJdpPtr_)
224 {
225 // If pointer is not set, set it to a zero field
226 dJdpPtr_.reset
227 (
228 createZeroFieldPtr<scalar>
229 (
230 mesh_,
231 ("dJdp_"+type()),
232 dimensionSet(0, 3, -2, 0, 0, 0, 0)
233 )
234 );
235 }
236 return *dJdpPtr_;
237}
238
239
241{
242 if (!dJdTPtr_)
243 {
244 // If pointer is not set, set it to a zero field
245 dJdTPtr_.reset
246 (
247 createZeroFieldPtr<scalar>
248 (
249 mesh_,
250 ("dJdT_"+type()),
251 dimensionSet(0, 3, -2, 0, 0, 0, 0)
252 )
253 );
254 }
255 return *dJdTPtr_;
256}
257
258
260{
261 if (!dJdTMvar1Ptr_)
262 {
263 // If pointer is not set, set it to a zero field
264 dJdTMvar1Ptr_.reset
265 (
266 createZeroFieldPtr<scalar>
267 (
268 mesh_,
269 ("dJdTMvar1_"+type()),
270 dimensionSet(0, 0, -2, 0, 0, 0, 0)
271 )
272 );
273 }
274 return *dJdTMvar1Ptr_;
275}
276
277
279{
280 if (!dJdTMvar2Ptr_)
281 {
282 // If pointer is not set, set it to a zero field
283 dJdTMvar2Ptr_.reset
284 (
285 createZeroFieldPtr<scalar>
286 (
287 mesh_,
288 ("dJdTMvar2_"+type()),
289 dimensionSet(0, 3, -2, 0, 0, 0, 0)
290 )
291 );
292 }
293 return *dJdTMvar2Ptr_;
294}
295
296
298(
299 const label patchI
300)
301{
302 if (!bdJdvPtr_)
303 {
304 bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
305 }
306 return bdJdvPtr_()[patchI];
307}
308
309
311(
312 const label patchI
313)
314{
315 if (!bdJdvnPtr_)
316 {
317 bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
318 }
319 return bdJdvnPtr_()[patchI];
320}
321
322
324(
325 const label patchI
326)
327{
328 if (!bdJdvtPtr_)
329 {
330 bdJdvtPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
331 }
332 return bdJdvtPtr_()[patchI];
333}
334
335
337(
338 const label patchI
339)
340{
341 if (!bdJdpPtr_)
342 {
343 bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
344 }
345 return bdJdpPtr_()[patchI];
346}
347
348
350(
351 const label patchI
352)
353{
354 if (!bdJdTPtr_)
355 {
356 bdJdTPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
357 }
358 return bdJdTPtr_()[patchI];
359}
360
361
363(
364 const label patchI
365)
366{
367 if (!bdJdTMvar1Ptr_)
368 {
369 bdJdTMvar1Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
370 }
371 return bdJdTMvar1Ptr_()[patchI];
372}
373
374
376(
377 const label patchI
378)
379{
380 if (!bdJdTMvar2Ptr_)
381 {
382 bdJdTMvar2Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
383 }
384 return bdJdTMvar2Ptr_()[patchI];
385}
386
387
389(
390 const label patchI
391)
392{
393 if (!bdJdnutPtr_)
394 {
395 bdJdnutPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
396 }
397 return bdJdnutPtr_()[patchI];
398}
399
400
402(
403 const label patchI
404)
405{
406 if (!bdJdGradUPtr_)
407 {
408 bdJdGradUPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
409 }
410 return bdJdGradUPtr_()[patchI];
411}
412
413
415{
416 if (!bdJdvPtr_)
417 {
418 bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
419 }
420 return bdJdvPtr_();
421}
422
423
425{
426 if (!bdJdvnPtr_)
427 {
428 bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
429 }
430 return bdJdvnPtr_();
431}
432
433
435{
436 if (!bdJdvtPtr_)
437 {
438 bdJdvtPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
439 }
440 return bdJdvtPtr_();
441}
442
443
445{
446 if (!bdJdpPtr_)
447 {
448 bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
449 }
450 return bdJdpPtr_();
451}
452
453
455{
456 if (!bdJdTPtr_)
457 {
458 bdJdTPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
459 }
460 return bdJdTPtr_();
461}
462
463
465{
466 if (!bdJdTMvar1Ptr_)
467 {
468 bdJdTMvar1Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
469 }
470 return bdJdTMvar1Ptr_();
471}
472
473
475{
476 if (!bdJdTMvar2Ptr_)
477 {
478 bdJdTMvar2Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
479 }
480 return bdJdTMvar2Ptr_();
481}
482
483
485{
486 if (!bdJdnutPtr_)
487 {
488 bdJdnutPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
489 }
490 return bdJdnutPtr_();
491}
492
493
495{
496 if (!bdJdGradUPtr_)
497 {
498 bdJdGradUPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
499 }
500 return *bdJdGradUPtr_;
501}
502
503
505{
506 // Objective function value
507 J();
508
509 // Update mean values here since they might be used in the
510 // subsequent functions
512
513 // volFields
514 update_dJdv();
515 update_dJdp();
516 update_dJdT();
519 update_dJdb();
522
523 // boundaryFields
539
540 // Divide everything with normalization factor
542}
543
544
546{
547 if (!nullified_)
548 {
549 if (hasdJdv())
550 {
551 dJdvPtr_() == dimensionedVector(dJdvPtr_().dimensions(), Zero);
552 }
553 if (hasdJdp())
554 {
555 dJdpPtr_() == dimensionedScalar(dJdpPtr_().dimensions(), Zero);
556 }
557 if (hasdJdT())
558 {
559 dJdTPtr_() == dimensionedScalar(dJdTPtr_().dimensions(), Zero);
560 }
561 if (hasdJdTMVar1())
562 {
563 dJdTMvar1Ptr_() ==
564 dimensionedScalar(dJdTMvar1Ptr_().dimensions(), Zero);
565 }
566 if (hasdJdTMVar2())
567 {
568 dJdTMvar2Ptr_() ==
569 dimensionedScalar(dJdTMvar2Ptr_().dimensions(), Zero);
570 }
571 if (hasBoundarydJdv())
572 {
574 }
575 if (hasBoundarydJdvn())
576 {
577 bdJdvnPtr_() == scalar(0);
578 }
579 if (hasBoundarydJdvt())
580 {
582 }
583 if (hasBoundarydJdp())
584 {
586 }
587 if (hasBoundarydJdT())
588 {
589 bdJdTPtr_() == scalar(0);
590 }
592 {
593 bdJdTMvar1Ptr_() == scalar(0);
594 }
596 {
597 bdJdTMvar2Ptr_() == scalar(0);
598 }
599 if (hasBoundarydJdnut())
600 {
601 bdJdnutPtr_() == scalar(0);
602 }
604 {
606 }
607
608 // Nullify geometric fields and sets nullified_ to true
610 }
611}
612
613
614bool objectiveIncompressible::write(const bool valid) const
615{
616 return objective::write(valid);
617}
618
619
620// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
621
622} // End namespace Foam
623
624// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Generic GeometricBoundaryField class.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
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
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for primal incompressible solvers.
void computeMeanFields()
Compute mean fields on the fly.
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.
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
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.
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 fvMesh & mesh_
Definition: objective.H:67
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:609
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:86
scalar weight_
Objective weight.
Definition: objective.H:83
bool computeMeanFields_
Definition: objective.H:72
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition: objective.C:336
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:94
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
Definition: objective.H:351
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
A non-counting (dummy) refCount.
Definition: refCount.H:59