adjointSimple.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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 FOSS GP
10  Copyright (C) 2019-2020 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 "adjointSimple.H"
31 #include "findRefCell.H"
32 #include "constrainHbyA.H"
33 #include "adjustPhi.H"
34 #include "fvOptions.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(adjointSimple, 0);
43  (
44  incompressibleAdjointSolver,
45  adjointSimple,
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 {
55  vars_.reset
56  (
58  (
59  mesh_,
63  )
64  );
65  return getAdjointVars();
66 }
67 
68 
70 {
71  if (adjointVars_.useSolverNameForFields())
72  {
74  << "useSolverNameForFields is set to true for adjointSolver "
75  << solverName() << nl << tab
76  << "Appending variable names with the solver name" << nl << tab
77  << "Please adjust the necessary entries in fvSchemes and fvSolution"
78  << nl << endl;
79  }
80 }
81 
82 
84 {
85  const surfaceScalarField& phia = adjointVars_.phiaInst();
87 
88  scalar sumLocalContErr = mesh_.time().deltaTValue()*
89  mag(contErr)().weightedAverage(mesh_.V()).value();
90 
91  scalar globalContErr = mesh_.time().deltaTValue()*
92  contErr.weightedAverage(mesh_.V()).value();
93  cumulativeContErr_ += globalContErr;
94 
95  Info<< "time step continuity errors : sum local = " << sumLocalContErr
96  << ", global = " << globalContErr
97  << ", cumulative = " << cumulativeContErr_
98  << endl;
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
104 Foam::adjointSimple::adjointSimple
105 (
106  fvMesh& mesh,
107  const word& managerType,
108  const dictionary& dict,
109  const word& primalSolverName
110 )
111 :
112  incompressibleAdjointSolver(mesh, managerType, dict, primalSolverName),
113  solverControl_(SIMPLEControl::New(mesh, managerType, *this)),
114  adjointVars_(allocateVars()),
115  cumulativeContErr_(Zero),
116  adjointSensitivity_(nullptr)
117 {
118  ATCModel_.reset
119  (
121  (
122  mesh,
123  primalVars_,
124  adjointVars_,
125  dict.subDict("ATCModel")
126  ).ptr()
127  );
128 
129  addExtraSchemes();
130  setRefCell
131  (
133  solverControl_().dict(),
136  );
137 
139  {
140  const IOdictionary& optDict =
141  mesh.lookupObject<IOdictionary>("optimisationDict");
142 
143  adjointSensitivity_.reset
144  (
146  (
147  mesh,
148  optDict.subDict("optimisation").subDict("sensitivities"),
149  primalVars_,
150  adjointVars_,
152  ).ptr()
153  );
154  }
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
160 bool Foam::adjointSimple::readDict(const dictionary& dict)
161 {
163  {
164  if (adjointSensitivity_)
165  {
166  const IOdictionary& optDict =
167  mesh_.lookupObject<IOdictionary>("optimisationDict");
168 
169  adjointSensitivity_().readDict
170  (
171  optDict.subDict("optimisation").subDict("sensitivities")
172  );
173  }
174 
175  return true;
176  }
177 
178  return false;
179 }
180 
181 
183 {
184  preIter();
185  mainIter();
186  postIter();
187 }
188 
189 
191 {
192  Info<< "Time = " << mesh_.time().timeName() << "\n" << endl;
193 }
194 
195 
197 {
198  // Grab primal references
199  const surfaceScalarField& phi = primalVars_.phi();
200  // Grab adjoint references
201  volScalarField& pa = adjointVars_.paInst();
202  volVectorField& Ua = adjointVars_.UaInst();
203  surfaceScalarField& phia = adjointVars_.phiaInst();
205  adjointVars_.adjointTurbulence();
206  const label& paRefCell = solverControl_().pRefCell();
207  const scalar& paRefValue = solverControl_().pRefValue();
208  fv::options& fvOptions(fv::options::New(this->mesh_));
209 
210  // Momentum predictor
211  //~~~~~~~~~~~~~~~~~~~
212 
213  tmp<fvVectorMatrix> tUaEqn
214  (
215  fvm::div(-phi, Ua)
216  + adjointTurbulence->divDevReff(Ua)
217  + adjointTurbulence->adjointMeanFlowSource()
218  ==
219  fvOptions(Ua)
220  );
221  fvVectorMatrix& UaEqn = tUaEqn.ref();
222 
223  // Add sources from boundary conditions
225 
226  // Add sources from volume-based objectives
227  objectiveManagerPtr_().addUaEqnSource(UaEqn);
228 
229  // Add ATC term
230  ATCModel_->addATC(UaEqn);
231 
232  // Additional source terms (e.g. energy equation)
233  addMomentumSource(UaEqn);
234 
235  UaEqn.relax();
236 
237  fvOptions.constrain(UaEqn);
238 
239  if (solverControl_().momentumPredictor())
240  {
241  Foam::solve(UaEqn == -fvc::grad(pa));
242 
243  fvOptions.correct(Ua);
244  }
245 
246  // Pressure Eq
247  //~~~~~~~~~~~~
248  {
249  volScalarField rAUa(1.0/UaEqn.A());
250  // 190402: Vag: to be updated.
251  // Probably a constrainHabyA by class is needed?
252  volVectorField HabyA(constrainHbyA(rAUa*UaEqn.H(), Ua, pa));
253  surfaceScalarField phiaHbyA("phiaHbyA", fvc::flux(HabyA));
254  adjustPhi(phiaHbyA, Ua, pa);
255 
256  tmp<volScalarField> rAtUa(rAUa);
257 
258  if (solverControl_().consistent())
259  {
260  rAtUa = 1.0/(1.0/rAUa - UaEqn.H1());
261  phiaHbyA +=
262  fvc::interpolate(rAtUa() - rAUa)*fvc::snGrad(pa)*mesh_.magSf();
263  HabyA -= (rAUa - rAtUa())*fvc::grad(pa);
264  }
265 
266  tUaEqn.clear();
267 
268  // Update the pressure BCs to ensure flux consistency
269  // constrainPressure(p, U, phiHbyA, rAtU(), MRF_);
270 
271  // Non-orthogonal pressure corrector loop
272  while (solverControl_().correctNonOrthogonal())
273  {
275  (
276  fvm::laplacian(rAtUa(), pa) == fvc::div(phiaHbyA)
277  );
278 
280 
281  addPressureSource(paEqn);
282 
284  paEqn.setReference(paRefCell, paRefValue);
285 
286  paEqn.solve();
287 
288  if (solverControl_().finalNonOrthogonalIter())
289  {
290  phia = phiaHbyA - paEqn.flux();
291  }
292  }
293 
294  continuityErrors();
295 
296  // Explicitly relax pressure for adjoint momentum corrector
297  pa.relax();
298 
299  // Momentum corrector
300  Ua = HabyA - rAtUa()*fvc::grad(pa);
302  fvOptions.correct(Ua);
304  }
305 
306  adjointTurbulence->correct();
307 
308  if (solverControl_().printMaxMags())
309  {
310  dimensionedScalar maxUa = gMax(mag(Ua)());
311  dimensionedScalar maxpa = gMax(mag(pa)());
312  Info<< "Max mag of adjoint velocity = " << maxUa.value() << endl;
313  Info<< "Max mag of adjoint pressure = " << maxpa.value() << endl;
314  }
315 }
316 
317 
319 {
320  solverControl_().write();
321 
322  // Average fields if necessary
323  adjointVars_.computeMeanFields();
324 
325  // Print execution time
326  mesh_.time().printExecutionTime(Info);
327 }
328 
329 
331 {
332  if (active_)
333  {
334  preLoop();
335  while (solverControl_().loop())
336  {
337  solveIter();
338  }
339  }
340 }
341 
342 
344 {
345  return solverControl_().loop();
346 }
347 
348 
350 {
351  // Reset mean fields before solving
352  adjointVars_.resetMeanFields();
353 }
354 
355 
357 {
358  if (computeSensitivities_)
359  {
360  adjointSensitivity_->accumulateIntegrand(scalar(1));
361  const scalarField& sens = adjointSensitivity_->calculateSensitivities();
362  if (!sensitivities_)
363  {
364  sensitivities_.reset(new scalarField(sens.size(), Zero));
365  }
366  *sensitivities_ = sens;
367  }
368  else
369  {
370  sensitivities_.reset(new scalarField());
371  }
372 }
373 
374 
376 {
377  if (!sensitivities_)
378  {
379  computeObjectiveSensitivities();
380  }
381 
382  return sensitivities_();
383 }
384 
385 
387 {
388  if (computeSensitivities_)
389  {
390  adjointSensitivity_->clearSensitivities();
392  }
393 }
394 
395 
397 {
398  if (!adjointSensitivity_.valid())
399  {
401  << "Sensitivity object not allocated" << nl
402  << "Turn computeSensitivities on in "
403  << solverName_
404  << nl << nl
406  }
407 
408  return adjointSensitivity_();
409 }
410 
411 
413 {
414  // Does nothing
415 }
416 
417 
419 {
420  // Does nothing
421 }
422 
423 
425 {
427 
428  // Update objective function related quantities
429  objectiveManagerPtr_->updateAndWrite();
430 }
431 
432 
434 {
435  os.writeEntry("averageIter", solverControl_().averageIter());
436 
437  return true;
438 }
439 
440 
441 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::adjointSimple::continuityErrors
void continuityErrors()
Definition: adjointSimple.C:76
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Foam::fvMatrix::H
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Definition: fvMatrix.C:1416
Foam::constrainHbyA
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:28
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
sumLocalContErr
scalar sumLocalContErr
Definition: compressibleContinuityErrs.H:28
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::fvMatrix::boundaryManipulate
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Definition: fvMatrix.C:1341
contErr
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Definition: adjustPhi.C:30
Foam::incompressibleAdjointSolver::ATCModel_
autoPtr< ATCModel > ATCModel_
Definition: incompressibleAdjointSolver.H:73
momentumPredictor
const bool momentumPredictor
Definition: readFluidMultiRegionSIMPLEControls.H:6
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Foam::adjointSolver::computeSensitivities_
bool computeSensitivities_
Definition: adjointSolver.H:80
Foam::adjointSimple::clearSensitivities
virtual void clearSensitivities()
Definition: adjointSimple.C:379
Foam::incompressibleAdjointSolver
Base class for incompressibleAdjoint solvers.
Definition: incompressibleAdjointSolver.H:47
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Definition: fvOptionListTemplates.C:348
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:280
Foam::adjointSimple::solveIter
virtual void solveIter()
Definition: adjointSimple.C:175
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::fv::optionList::constrain
void constrain(fvMatrix< Type > &eqn)
Definition: fvOptionListTemplates.C:307
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Definition: fvOptions.C:96
findRefCell.H
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
Foam::adjointSimple::computeObjectiveSensitivities
virtual void computeObjectiveSensitivities()
Definition: adjointSimple.C:349
Foam::adjointSimple::addMomentumSource
virtual void addMomentumSource(fvVectorMatrix &matrix)
Definition: adjointSimple.C:405
Foam::solver::mesh_
fvMesh & mesh_
Definition: solver.H:54
Foam::adjointSimple::adjointVars_
incompressibleAdjointVars & adjointVars_
Definition: adjointSimple.H:84
Foam::adjointSimple::preLoop
virtual void preLoop()
Definition: adjointSimple.C:342
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::dimensioned::value
const Type & value() const
Definition: dimensionedType.C:427
Foam::fvMatrix::setReference
void setReference(const label celli, const Type &value, const bool forceReference=false)
Definition: fvMatrix.C:1085
adjustPhi.H
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:46
Foam::adjointSolver::objectiveManagerPtr_
autoPtr< objectiveManager > objectiveManagerPtr_
Definition: adjointSolver.H:74
Foam::adjointSimple::preIter
virtual void preIter()
Definition: adjointSimple.C:183
Foam::adjointSimple::adjointSensitivity_
autoPtr< incompressible::adjointSensitivity > adjointSensitivity_
Definition: adjointSimple.H:90
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Foam::incompressibleAdjointMeanFlowVars::paInst
const volScalarField & paInst() const
Definition: incompressibleAdjointMeanFlowVars.C:218
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:220
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:38
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::adjointSimple::addPressureSource
virtual void addPressureSource(fvScalarMatrix &matrix)
Definition: adjointSimple.C:411
Foam::Info
messageStream Info
Foam::adjointSimple::getSensitivityBase
virtual sensitivity & getSensitivityBase()
Definition: adjointSimple.C:389
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Foam::incompressible::adjointSensitivity::New
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, incompressibleVars &primalVars, incompressibleAdjointVars &adjointVars, objectiveManager &objectiveManager)
Definition: adjointSensitivityIncompressible.C:65
Foam::fvMatrix::H1
tmp< volScalarField > H1() const
Definition: fvMatrix.C:1478
Foam::fvMatrix::relax
void relax(const scalar alpha)
Definition: fvMatrix.C:1176
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:453
Foam::setRefCell
bool setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
Definition: findRefCell.C:27
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:427
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:51
paEqn
fvScalarMatrix paEqn(fvm::d2dt2(pa) - sqr(c0) *fvc::laplacian(pa))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::incompressibleAdjointSolver::updatePrimalBasedQuantities
virtual void updatePrimalBasedQuantities()
Definition: incompressibleAdjointSolver.C:153
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::adjointSolver::clearSensitivities
virtual void clearSensitivities()
Definition: adjointSolver.C:153
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:36
Foam::fvMatrix::A
tmp< volScalarField > A() const
Definition: fvMatrix.C:1387
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam
Definition: atmBoundaryLayer.C:26
Foam::fvc::interpolate
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Foam::adjointSimple::solve
virtual void solve()
Definition: adjointSimple.C:323
Foam::SIMPLEControl::New
static autoPtr< SIMPLEControl > New(fvMesh &mesh, const word &managerType, const solver &solver)
Definition: SIMPLEControl.C:57
globalContErr
scalar globalContErr
Definition: compressibleContinuityErrs.H:31
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::adjointSimple::mainIter
virtual void mainIter()
Definition: adjointSimple.C:189
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Definition: GeometricField.C:933
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Definition: DimensionedFieldReuseFunctions.H:100
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::solver::mesh
const fvMesh & mesh() const
Definition: solver.C:89
Foam::adjointSimple::readDict
virtual bool readDict(const dictionary &dict)
Definition: adjointSimple.C:153
Foam::tab
constexpr char tab
Definition: Ostream.H:423
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::solver::vars_
autoPtr< variablesSet > vars_
Definition: solver.H:76
Foam::DimensionedField::weightedAverage
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Definition: DimensionedField.C:454
Foam::adjointSimple::updatePrimalBasedQuantities
virtual void updatePrimalBasedQuantities()
Definition: adjointSimple.C:417
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:50
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Definition: GeometricField.C:776
Foam::adjointSimple::postIter
virtual void postIter()
Definition: adjointSimple.C:311
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
adjointSimple.H
Foam::incompressibleAdjointSolver::primalVars_
incompressibleVars & primalVars_
Definition: incompressibleAdjointSolver.H:70
Foam::GeometricField::relax
void relax(const scalar alpha)
Definition: GeometricField.C:965
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:64
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Definition: Ostream.H:232
Foam::adjointSimple::writeData
virtual bool writeData(Ostream &os) const
Definition: adjointSimple.C:426
Foam::incompressibleAdjointSolver::getAdjointVars
virtual const incompressibleAdjointVars & getAdjointVars() const
Definition: incompressibleAdjointSolver.C:122
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Definition: fvcFlux.C:27
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::adjointSimple::getObjectiveSensitivities
virtual const scalarField & getObjectiveSensitivities()
Definition: adjointSimple.C:368
Foam::sensitivity
Abstract base class for adjoint sensitivities.
Definition: sensitivity.H:59
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Foam::incompressibleAdjointSolver::readDict
virtual bool readDict(const dictionary &dict)
Definition: incompressibleAdjointSolver.C:97
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::adjointSimple::addExtraSchemes
void addExtraSchemes()
Definition: adjointSimple.C:62
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
constrainHbyA.H
Foam::fvMatrix::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Definition: fvMatrix.C:1527
Foam::solver::dict
virtual const dictionary & dict() const
Definition: solver.C:106
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:585
Foam::adjointSimple::loop
virtual bool loop()
Definition: adjointSimple.C:336
Foam::adjointSimple::solverControl_
autoPtr< SIMPLEControl > solverControl_
Definition: adjointSimple.H:79
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:36
Foam::ATCModel::New
static autoPtr< ATCModel > New(const fvMesh &mesh, const incompressibleVars &primalVars, const incompressibleAdjointVars &adjointVars, const dictionary &dict)
Definition: ATCModel.C:122
Foam::SIMPLEControl
SIMPLE control class to supply convergence information/checks for the SIMPLE loop.
Definition: SIMPLEControl.H:46
Foam::adjointSimple::allocateVars
incompressibleAdjointVars & allocateVars()
Definition: adjointSimple.C:46
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:37