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 -------------------------------------------------------------------------------
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 
33 #include "createZeroField.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(objectiveIncompressible, 0);
44 defineRunTimeSelectionTable(objectiveIncompressible, dictionary);
46 (
47  objective,
49  objective
50 );
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 objectiveIncompressible::objectiveIncompressible
55 (
56  const fvMesh& mesh,
57  const dictionary& dict,
58  const word& adjointSolverName,
59  const word& primalSolverName
60 )
61 :
62  objective(mesh, dict, adjointSolverName, primalSolverName),
63 
64  vars_
65  (
66  mesh.lookupObject<incompressiblePrimalSolver>(primalSolverName).
67  getIncoVars()
68  ),
69 
70  // Initialize pointers to nullptr.
71  // Not all of them are required for each objective function.
72  // Each child should allocate whatever is needed.
73 
74  // Field adjoint Eqs
75  dJdvPtr_(nullptr),
76  dJdpPtr_(nullptr),
77  dJdTPtr_(nullptr),
78  dJdTMvar1Ptr_(nullptr),
79  dJdTMvar2Ptr_(nullptr),
80 
81  // Adjoint boundary conditions
82  bdJdvPtr_(nullptr),
83  bdJdvnPtr_(nullptr),
84  bdJdvtPtr_(nullptr),
85  bdJdpPtr_(nullptr),
86  bdJdTPtr_(nullptr),
87  bdJdTMvar1Ptr_(nullptr),
88  bdJdTMvar2Ptr_(nullptr),
89  bdJdnutPtr_(nullptr),
90  bdJdGradUPtr_(nullptr)
91 {
92  weight_ = dict.get<scalar>("weight");
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
98 
100 (
101  const fvMesh& mesh,
102  const dictionary& dict,
103  const word& adjointSolverName,
104  const word& primalSolverName
105 )
106 {
107  const word modelType(dict.get<word>("type"));
108 
109  Info<< "Creating objective function : " << dict.dictName()
110  << " of type " << modelType << endl;
111 
112  auto* ctorPtr = dictionaryConstructorTable(modelType);
113 
114  if (!ctorPtr)
115  {
117  (
118  dict,
119  "objectiveIncompressible",
120  modelType,
121  *dictionaryConstructorTablePtr_
122  ) << exit(FatalIOError);
123  }
124 
125  return autoPtr<objectiveIncompressible>
126  (
127  ctorPtr(mesh, dict, adjointSolverName, primalSolverName)
128  );
129 }
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 {
136  if (normalize_ && normFactor_)
137  {
138  const scalar oneOverNorm(1./normFactor_());
139 
140  if (hasdJdv())
141  {
142  dJdvPtr_().primitiveFieldRef() *= oneOverNorm;
143  }
144  if (hasdJdp())
145  {
146  dJdpPtr_().primitiveFieldRef() *= oneOverNorm;
147  }
148  if (hasdJdT())
149  {
150  dJdTPtr_().primitiveFieldRef() *= oneOverNorm;
151  }
152  if (hasdJdTMVar1())
153  {
154  dJdTMvar1Ptr_().primitiveFieldRef() *= oneOverNorm;
155  }
156  if (hasdJdTMVar2())
157  {
158  dJdTMvar2Ptr_().primitiveFieldRef() *= oneOverNorm;
159  }
160  if (hasBoundarydJdv())
161  {
162  bdJdvPtr_() *= oneOverNorm;
163  }
164  if (hasBoundarydJdvn())
165  {
166  bdJdvnPtr_() *= oneOverNorm;
167  }
168  if (hasBoundarydJdvt())
169  {
170  bdJdvtPtr_() *= oneOverNorm;
171  }
172  if (hasBoundarydJdp())
173  {
174  bdJdpPtr_() *= oneOverNorm;
175  }
176  if (hasBoundarydJdT())
177  {
178  bdJdTPtr_() *= oneOverNorm;
179  }
180  if (hasBoundarydJdTMVar1())
181  {
182  bdJdTMvar1Ptr_() *= oneOverNorm;
183  }
184  if (hasBoundarydJdTMVar2())
185  {
186  bdJdTMvar2Ptr_() *= oneOverNorm;
187  }
188  if (hasBoundarydJdnut())
189  {
190  bdJdnutPtr_() *= oneOverNorm;
191  }
192  if (hasBoundarydJdGradU())
193  {
194  bdJdGradUPtr_() *= oneOverNorm;
195  }
197  // Normalize geometric fields
199  }
200 }
201 
202 
204 {
205  if (!dJdvPtr_)
206  {
207  // If pointer is not set, set it to a zero field
208  dJdvPtr_.reset
209  (
210  createZeroFieldPtr<vector>
211  (
212  mesh_,
213  ("dJdv_"+type()),
215  )
216  );
217  }
218  return *dJdvPtr_;
219 }
220 
221 
223 {
224  if (!dJdpPtr_)
225  {
226  // If pointer is not set, set it to a zero field
227  dJdpPtr_.reset
228  (
229  createZeroFieldPtr<scalar>
230  (
231  mesh_,
232  ("dJdp_"+type()),
233  dimensionSet(0, 3, -2, 0, 0, 0, 0)
234  )
235  );
236  }
237  return *dJdpPtr_;
238 }
239 
240 
242 {
243  if (!dJdTPtr_)
244  {
245  // If pointer is not set, set it to a zero field
246  dJdTPtr_.reset
247  (
248  createZeroFieldPtr<scalar>
249  (
250  mesh_,
251  ("dJdT_"+type()),
252  dimensionSet(0, 3, -2, 0, 0, 0, 0)
253  )
254  );
255  }
256  return *dJdTPtr_;
257 }
258 
259 
261 {
262  if (!dJdTMvar1Ptr_)
263  {
264  // If pointer is not set, set it to a zero field
265  dJdTMvar1Ptr_.reset
266  (
267  createZeroFieldPtr<scalar>
268  (
269  mesh_,
270  ("dJdTMvar1_"+type()),
271  dimensionSet(0, 0, -2, 0, 0, 0, 0)
272  )
273  );
274  }
275  return *dJdTMvar1Ptr_;
276 }
277 
278 
280 {
281  if (!dJdTMvar2Ptr_)
282  {
283  // If pointer is not set, set it to a zero field
284  dJdTMvar2Ptr_.reset
285  (
286  createZeroFieldPtr<scalar>
287  (
288  mesh_,
289  ("dJdTMvar2_"+type()),
290  dimensionSet(0, 3, -2, 0, 0, 0, 0)
291  )
292  );
293  }
294  return *dJdTMvar2Ptr_;
295 }
296 
297 
299 (
300  const label patchI
301 )
302 {
303  if (!bdJdvPtr_)
304  {
305  bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
306  }
307  return bdJdvPtr_()[patchI];
308 }
309 
310 
312 (
313  const label patchI
314 )
315 {
316  if (!bdJdvnPtr_)
317  {
318  bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
319  }
320  return bdJdvnPtr_()[patchI];
321 }
322 
323 
325 (
326  const label patchI
327 )
328 {
329  if (!bdJdvtPtr_)
330  {
331  bdJdvtPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
332  }
333  return bdJdvtPtr_()[patchI];
334 }
335 
336 
338 (
339  const label patchI
340 )
341 {
342  if (!bdJdpPtr_)
343  {
344  bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
345  }
346  return bdJdpPtr_()[patchI];
347 }
348 
349 
351 (
352  const label patchI
353 )
354 {
355  if (!bdJdTPtr_)
356  {
357  bdJdTPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
358  }
359  return bdJdTPtr_()[patchI];
360 }
361 
362 
364 (
365  const label patchI
366 )
367 {
368  if (!bdJdTMvar1Ptr_)
369  {
370  bdJdTMvar1Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
371  }
372  return bdJdTMvar1Ptr_()[patchI];
373 }
374 
375 
377 (
378  const label patchI
379 )
380 {
381  if (!bdJdTMvar2Ptr_)
382  {
383  bdJdTMvar2Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
384  }
385  return bdJdTMvar2Ptr_()[patchI];
386 }
387 
388 
390 (
391  const label patchI
392 )
393 {
394  if (!bdJdnutPtr_)
395  {
396  bdJdnutPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
397  }
398  return bdJdnutPtr_()[patchI];
399 }
400 
401 
403 (
404  const label patchI
405 )
406 {
407  if (!bdJdGradUPtr_)
408  {
409  bdJdGradUPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
410  }
411  return bdJdGradUPtr_()[patchI];
412 }
413 
414 
416 {
417  if (!bdJdvPtr_)
418  {
419  bdJdvPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
420  }
421  return bdJdvPtr_();
422 }
423 
424 
426 {
427  if (!bdJdvnPtr_)
428  {
429  bdJdvnPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
430  }
431  return bdJdvnPtr_();
432 }
433 
434 
436 {
437  if (!bdJdvtPtr_)
438  {
439  bdJdvtPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
440  }
441  return bdJdvtPtr_();
442 }
443 
444 
446 {
447  if (!bdJdpPtr_)
448  {
449  bdJdpPtr_.reset(createZeroBoundaryPtr<vector>(mesh_));
450  }
451  return bdJdpPtr_();
452 }
453 
454 
456 {
457  if (!bdJdTPtr_)
458  {
459  bdJdTPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
460  }
461  return bdJdTPtr_();
462 }
463 
464 
466 {
467  if (!bdJdTMvar1Ptr_)
468  {
469  bdJdTMvar1Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
470  }
471  return bdJdTMvar1Ptr_();
472 }
473 
474 
476 {
477  if (!bdJdTMvar2Ptr_)
478  {
479  bdJdTMvar2Ptr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
480  }
481  return bdJdTMvar2Ptr_();
482 }
483 
484 
486 {
487  if (!bdJdnutPtr_)
488  {
489  bdJdnutPtr_.reset(createZeroBoundaryPtr<scalar>(mesh_));
490  }
491  return bdJdnutPtr_();
492 }
493 
494 
496 {
497  if (!bdJdGradUPtr_)
498  {
499  bdJdGradUPtr_.reset(createZeroBoundaryPtr<tensor>(mesh_));
500  }
501  return *bdJdGradUPtr_;
502 }
503 
504 
506 {
507  // Objective function value
508  J();
509 
510  // Update mean values here since they might be used in the
511  // subsequent functions
513 
514  // volFields
515  update_dJdv();
516  update_dJdp();
517  update_dJdT();
520  update_dJdb();
523 
524  // boundaryFields
540 
541  // Divide everything with normalization factor
542  doNormalization();
543 }
544 
545 
547 {
548  if (!nullified_)
549  {
550  if (hasdJdv())
551  {
552  dJdvPtr_() == dimensionedVector(dJdvPtr_().dimensions(), Zero);
553  }
554  if (hasdJdp())
555  {
556  dJdpPtr_() == dimensionedScalar(dJdpPtr_().dimensions(), Zero);
557  }
558  if (hasdJdT())
559  {
560  dJdTPtr_() == dimensionedScalar(dJdTPtr_().dimensions(), Zero);
561  }
562  if (hasdJdTMVar1())
563  {
564  dJdTMvar1Ptr_() ==
565  dimensionedScalar(dJdTMvar1Ptr_().dimensions(), Zero);
566  }
567  if (hasdJdTMVar2())
568  {
569  dJdTMvar2Ptr_() ==
570  dimensionedScalar(dJdTMvar2Ptr_().dimensions(), Zero);
571  }
572  if (hasBoundarydJdv())
573  {
574  bdJdvPtr_() == vector::zero;
575  }
576  if (hasBoundarydJdvn())
577  {
578  bdJdvnPtr_() == scalar(0);
579  }
580  if (hasBoundarydJdvt())
581  {
582  bdJdvtPtr_() == vector::zero;
583  }
584  if (hasBoundarydJdp())
585  {
586  bdJdpPtr_() == vector::zero;
587  }
588  if (hasBoundarydJdT())
589  {
590  bdJdTPtr_() == scalar(0);
591  }
592  if (hasBoundarydJdTMVar1())
593  {
594  bdJdTMvar1Ptr_() == scalar(0);
595  }
596  if (hasBoundarydJdTMVar2())
597  {
598  bdJdTMvar2Ptr_() == scalar(0);
599  }
600  if (hasBoundarydJdnut())
601  {
602  bdJdnutPtr_() == scalar(0);
603  }
604  if (hasBoundarydJdGradU())
605  {
606  bdJdGradUPtr_() == tensor::zero;
607  }
609  // Nullify geometric fields and sets nullified_ to true
611  }
612 }
613 
614 
616 {
617  // Figure out the availability of the adjoint turbulence model variables
618  // through their primal counterparts, since the contructor of the adjoint
619  // solver has not been completed yet
620  const incompressiblePrimalSolver& primSolver =
621  mesh_.lookupObject<incompressiblePrimalSolver>(primalSolverName_);
622  const autoPtr<incompressible::RASModelVariables>& rasVars =
623  primSolver.getIncoVars().RASModelVariables();
624 
625  if (rasVars().hasTMVar1())
626  {
627  const dimensionSet primalVarDims = rasVars->TMVar1Inst().dimensions();
628  const dimensionSet sourceDims = dimArea/pow3(dimTime)/primalVarDims;
629  dJdTMvar1Ptr_.reset
630  (
631  createZeroFieldPtr<scalar>
632  (
633  mesh_,
634  "ATMSource1",
635  sourceDims
636  )
637  );
638  }
639  if (rasVars().hasTMVar2())
640  {
641  const dimensionSet primalVarDims = rasVars->TMVar2Inst().dimensions();
642  const dimensionSet sourceDims = dimArea/pow3(dimTime)/primalVarDims;
643  dJdTMvar2Ptr_.reset
644  (
645  createZeroFieldPtr<scalar>
646  (
647  mesh_,
648  "ATMSource2",
649  sourceDims
650  )
651  );
652  }
653 }
654 
655 
657 (
658  const labelList& zoneIDs
659 ) const
660 {
661  label nCells(0);
662  for (const label zI : zoneIDs)
663  {
664  nCells += mesh_.cellZones()[zI].size();
665  }
666  reduce(nCells, sumOp<label>());
667  if (!nCells)
668  {
670  << "Provided cellZones include no cells"
671  << exit(FatalError);
672  }
673 }
674 
675 
677 (
678  autoPtr<volScalarField>& dJdTMvarPtr,
679  tmp<volScalarField>
680  (incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const,
681  const volScalarField& JacobianMultiplier,
682  const labelList& zones
683 )
684 {
685  if (dJdTMvarPtr)
686  {
687  // nut Jacobians are currently computed in the adjoint turbulence
688  // models, though they would be better placed within the primal
689  // turbulence model.
690  // Safeguard the computation until the adjoint solver is complete
691  if (mesh_.foundObject<incompressibleAdjointSolver>(adjointSolverName_))
692  {
693  const incompressibleAdjointSolver& adjSolver =
694  mesh_.lookupObject<incompressibleAdjointSolver>
696  const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointRAS =
697  adjSolver.getAdjointVars().adjointTurbulence();
698 
699  tmp<volScalarField> tnutJacobian = (adjointRAS->*JacobianFunc)();
700  const volScalarField& nutJacobian = tnutJacobian();
701 
702  volScalarField& dJdTMvar = dJdTMvarPtr();
703 
704  for (const label zI : zones)
705  {
706  const cellZone& zoneI = mesh_.cellZones()[zI];
707  for (const label cellI : zoneI)
708  {
709  dJdTMvar[cellI] =
710  JacobianMultiplier[cellI]*nutJacobian[cellI];
711  }
712  }
713  }
714  else
715  {
717  << "Skipping the computation of nutJacobian until "
718  << "the adjoint solver is complete"
719  << endl;
720  }
721  }
722 }
723 
724 
725 bool objectiveIncompressible::write(const bool valid) const
726 {
727  return objective::write(valid);
728 }
729 
730 
731 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
732 
733 } // End namespace Foam
734 
735 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
virtual void update_meanValues()
Some objectives need to store some auxiliary values. If averaging is enabled, update these mean value...
dictionary dict
virtual void doNormalization()
Normalize all fields allocated by the objective.
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 labelIOList & zoneIDs
Definition: correctPhi.H:59
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const volScalarField & dJdp()
Contribution to field adjoint continuity eq.
fvPatchField< vector > fvPatchVectorField
volTensorField::Boundary boundaryTensorField
const word adjointSolverName_
Definition: objective.H:65
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:87
autoPtr< volScalarField > dJdpPtr_
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
autoPtr< boundaryScalarField > bdJdTMvar2Ptr_
Adjoint outlet turbulence model var 2.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
const boundaryScalarField & boundarydJdT()
Objective partial deriv wrt temperature for all patches.
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
virtual void update()
Update objective function derivatives.
autoPtr< volScalarField > dJdTPtr_
dimensionedSymmTensor sqr(const dimensionedVector &dv)
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:90
volVectorField::Boundary boundaryVectorField
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition: objective.C:329
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
fvPatchField< tensor > fvPatchTensorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
Macros for easy insertion into run-time selection tables.
const boundaryScalarField & boundarydJdTMvar2()
Objective partial deriv wrt turbulence model var 2 for all patches.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
const volScalarField & dJdT()
Contribution to field adjoint energy eq.
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
virtual bool write(const bool valid=true) const
Write objective function history.
Definition: objective.C:653
const boundaryScalarField & boundarydJdTMvar1()
Objective partial deriv wrt turbulence model var 1 for all patches.
dynamicFvMesh & mesh
const boundaryScalarField & boundarydJdvn()
Objective partial deriv wrt normal velocity for all patches.
bool computeMeanFields_
Definition: objective.H:68
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
const boundaryTensorField & boundarydJdGradU()
Objective partial deriv wrt gradU.
virtual void nullify()
Update objective function derivatives.
autoPtr< volVectorField > dJdvPtr_
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
const volScalarField & dJdTMvar1()
Contribution to field adjoint turbulence model variable 1.
const fvMesh & mesh_
Definition: objective.H:63
const boundaryVectorField & boundarydJdp()
Objective partial deriv wrt pressure (times normal) for all patches.
const volScalarField & dJdTMvar2()
Contribution to field adjoint turbulence model variable 2.
virtual scalar J()=0
Return the objective function value.
autoPtr< boundaryVectorField > bdJdvPtr_
virtual bool write(const bool valid=true) const
Write objective function history.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:55
const boundaryVectorField & boundarydJdvt()
Objective partial deriv wrt tangent velocity for all patches.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
defineTypeNameAndDebug(combustionModel, 0)
bool hasdJdv() const
Inline functions for checking whether pointers are set or not.
const word primalSolverName_
Definition: objective.H:66
autoPtr< boundaryScalarField > bdJdTPtr_
Adjoint outlet temperature.
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
void computeMeanFields()
Compute mean fields on the fly.
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
const incompressibleVars & vars_
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
Definition: objective.H:479
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:602
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
const boundaryScalarField & boundarydJdnut()
Objective partial deriv wrt nut for all patches.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar weight_
Objective weight.
Definition: objective.H:85
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
virtual void update_dJdv()
Update vol and boundary fields and derivatives.
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
virtual void update_dxdbDirectMultiplier()
Update d (x) / db multiplier. Surface and volume-based sensitivity term.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:654
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
autoPtr< boundaryScalarField > bdJdTMvar1Ptr_
Adjoint outlet turbulence model var 1.
Base class for primal incompressible solvers.
List< label > labelList
A List of labels.
Definition: List.H:62
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
virtual void update_boundarydJdb()
Update objective function derivative term.
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:615
volScalarField::Boundary boundaryScalarField
Abstract base class for objective functions in incompressible flows.
const boundaryVectorField & boundarydJdv()
Objective partial deriv wrt velocity for all patches.
const volVectorField & dJdv()
Contribution to field adjoint momentum eqs.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157