Loading...
Searching...
No Matches
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-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\*---------------------------------------------------------------------------*/
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
39namespace Foam
40{
43 (
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 const surfaceScalarField& phia = adjointVars_.phiaInst();
73
74 scalar sumLocalContErr = mesh_.time().deltaTValue()*
75 mag(contErr)().weightedAverage(mesh_.V()).value();
76
77 scalar globalContErr = mesh_.time().deltaTValue()*
78 contErr.weightedAverage(mesh_.V()).value();
79 cumulativeContErr_ += globalContErr;
80
81 Info<< "time step continuity errors : sum local = " << sumLocalContErr
82 << ", global = " << globalContErr
83 << ", cumulative = " << cumulativeContErr_
84 << endl;
85}
86
87
90 adjointSensitivity_->accumulateIntegrand(scalar(1));
91}
92
93
94// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95
96Foam::adjointSimple::adjointSimple
97(
98 fvMesh& mesh,
99 const word& managerType,
100 const dictionary& dict,
101 const word& primalSolverName,
102 const word& solverName
103)
104:
106 (
107 mesh,
108 managerType,
109 dict,
110 primalSolverName,
111 solverName
112 ),
113 solverControl_(SIMPLEControl::New(mesh, managerType, *this)),
114 adjointVars_(allocateVars()),
115 cumulativeContErr_(Zero)
116{
117 ATCModel_.reset
118 (
120 (
121 mesh,
124 dict.subDict("ATCModel")
125 ).ptr()
126 );
127
129 (
134 );
136}
137
138
139// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140
142{
143 solverControl_().incrementIter();
144 if (solverControl_().performIter())
145 {
147 mainIter();
148 postIter();
149 }
150}
151
154{
155 Info<< "Time = " << mesh_.time().timeName() << "\n" << endl;
156}
157
158
160{
161 addProfiling(adjointSimple, "adjointSimple::mainIter");
162 // Grab primal references
163 const surfaceScalarField& phi = primalVars_.phi();
164 // Grab adjoint references
165 volScalarField& pa = adjointVars_.paInst();
166 volVectorField& Ua = adjointVars_.UaInst();
167 surfaceScalarField& phia = adjointVars_.phiaInst();
169 adjointVars_.adjointTurbulence();
170 const label& paRefCell = solverControl_().pRefCell();
171 const scalar& paRefValue = solverControl_().pRefValue();
173
174 // Momentum predictor
175 //~~~~~~~~~~~~~~~~~~~
176
178 (
179 fvm::div(-phi, Ua)
180 + adjointTurbulence->divDevReff(Ua)
181 + adjointTurbulence->adjointMeanFlowSource()
182 ==
183 fvOptions(Ua)
184 );
185 fvVectorMatrix& UaEqn = tUaEqn.ref();
186
187 // Add sources from boundary conditions
189
190 // Add sources from volume-based objectives
191 objectiveManager_.addSource(UaEqn);
192
193 // Add ATC term
194 ATCModel_->addATC(UaEqn);
195
196 // Additional source terms (e.g. energy equation)
197 addMomentumSource(UaEqn);
198
199 UaEqn.relax();
200
201 fvOptions.constrain(UaEqn);
202
203 if (solverControl_().momentumPredictor())
204 {
205 Foam::solve(UaEqn == -fvc::grad(pa));
206
207 fvOptions.correct(Ua);
208 }
209
210 // Pressure Eq
211 //~~~~~~~~~~~~
212 {
213 volScalarField rAUa(1.0/UaEqn.A());
214 // 190402: Vag: to be updated.
215 // Probably a constrainHabyA by class is needed?
216 volVectorField HabyA(constrainHbyA(rAUa*UaEqn.H(), Ua, pa));
217 surfaceScalarField phiaHbyA("phiaHbyA", fvc::flux(HabyA));
218 adjustPhi(phiaHbyA, Ua, pa);
219
220 tmp<volScalarField> rAtUa(rAUa);
221
222 if (solverControl_().consistent())
223 {
224 rAtUa = 1.0/(1.0/rAUa - UaEqn.H1());
225 phiaHbyA +=
226 fvc::interpolate(rAtUa() - rAUa)*fvc::snGrad(pa)*mesh_.magSf();
227 HabyA -= (rAUa - rAtUa())*fvc::grad(pa);
228 }
229
230 tUaEqn.clear();
231
232 // Update the pressure BCs to ensure flux consistency
233 // constrainPressure(p, U, phiHbyA, rAtU(), MRF_);
234
235 // Non-orthogonal pressure corrector loop
236 while (solverControl_().correctNonOrthogonal())
237 {
239 (
240 fvm::laplacian(rAtUa(), pa) == fvc::div(phiaHbyA)
241 );
242
244
245 addPressureSource(paEqn);
246
248 paEqn.setReference(paRefCell, paRefValue);
249
250 paEqn.solve();
251
252 if (solverControl_().finalNonOrthogonalIter())
253 {
254 phia = phiaHbyA - paEqn.flux();
255 }
256 }
257
258 continuityErrors();
259
260 // Explicitly relax pressure for adjoint momentum corrector
261 pa.relax();
262
263 // Momentum corrector
264 Ua = HabyA - rAtUa()*fvc::grad(pa);
266 fvOptions.correct(Ua);
268 }
269
270 adjointTurbulence->correct();
271
272 if (solverControl_().printMaxMags())
273 {
274 dimensionedScalar maxUa = gMax(mag(Ua)());
275 dimensionedScalar maxpa = gMax(mag(pa)());
276 Info<< "Max mag (" << Ua.name() << ") = " << maxUa.value() << endl;
277 Info<< "Max mag (" << pa.name() << ") = " << maxpa.value() << endl;
278 }
279}
280
281
283{
284 solverControl_().write();
285
286 // Average fields if necessary
287 adjointVars_.computeMeanFields();
288
289 // Print execution time
290 mesh_.time().printExecutionTime(Info);
291}
292
293
295{
296 addProfiling(adjointSimple, "adjointSimple::solve");
297 if (active_)
298 {
299 preLoop();
300 while (solverControl_().loop())
301 {
303 }
304 postLoop();
305 }
306}
307
310{
311 return solverControl_().loop();
312}
313
314
317 // Reset initial and mean fields before solving
318 adjointVars_.restoreInitValues();
319 adjointVars_.resetMeanFields();
320}
321
324{
325 // Does nothing
326}
327
330{
331 // Does nothing
332}
333
334
336{
338
339 // Update objective function related quantities
340 objectiveManager_.updateAndWrite();
341}
342
343
345{
346 // Determine number of variables related to the adjoint turbulence model
348 adjointVars_.adjointTurbulence();
349 const wordList& turbVarNames =
350 adjointTurbulence().getAdjointTMVariablesBaseNames();
351 label nTurbVars(turbVarNames.size());
352 if (adjointTurbulence().includeDistance())
353 {
354 nTurbVars++;
355 }
356
357 // Determine names of fields to be added to the dictionary
358 wordList names(1 + nTurbVars);
359 label varID(0);
360 names[varID++] = adjointVars_.UaInst().name();
361 for (const word& turbName : turbVarNames)
362 {
363 names[varID++] = turbName;
364 }
365 if (adjointTurbulence().includeDistance())
366 {
367 names[varID++] =
368 word(useSolverNameForFields() ? "da" + solverName_ : "da");
369 }
370
371 // Add entries to dictionary
372 const word dictName("topOSource" + solverName_);
373 dictionary optionDict(dictName);
374 optionDict.add<word>("type", "topOSource");
375 optionDict.add<wordList>("names", names);
376 optionDict.add<word>("function", "linear");
377 optionDict.add<word>("interpolationField", "beta");
379 // Construct and append fvOption
381 fvOptions.push_back(fv::option::New(dictName, optionDict, mesh_));
382}
383
384
386{
387 os.writeEntry("averageIter", solverControl_().averageIter());
389}
390
391
392// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
fv::options & fvOptions
const scalar pRefValue
const label pRefCell
static autoPtr< ATCModel > New(const fvMesh &mesh, const incompressibleVars &primalVars, const incompressibleAdjointVars &adjointVars, const dictionary &dict)
Return a reference to the selected turbulence model.
Definition ATCModel.C:123
void relax(const scalar alpha)
Relax field (for steady-state solution).
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
SIMPLE control class to supply convergence information/checks for the SIMPLE loop.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Solution of the adjoint PDEs for incompressible, steady-state flows.
virtual void postIter()
Steps to be executed before each main SIMPLE iteration.
autoPtr< SIMPLEControl > solverControl_
Solver control.
scalar cumulativeContErr_
Cumulative continuity error.
virtual void preIter()
Steps to be executed before each main SIMPLE iteration.
virtual bool writeData(Ostream &os) const
Write average iteration.
virtual void addMomentumSource(fvVectorMatrix &matrix)
Source terms for the momentum equations.
incompressibleAdjointVars & allocateVars()
Allocate incompressibleAdjointVars and return reference to be used for convenience in the rest of the...
virtual void mainIter()
The main SIMPLE iter.
virtual void addTopOFvOptions() const
Add fvOptions for TopO.
virtual void updatePrimalBasedQuantities()
Update primal based quantities related to the objective functions.
virtual void preCalculateSensitivities()
Accumulate the sensitivities integrand before calculating them.
incompressibleAdjointVars & adjointVars_
Reference to incompressibleAdjointVars.
virtual void preLoop()
Functions to be called before loop.
void continuityErrors()
Compute continuity errors.
virtual bool loop()
Looper (advances iters, time step).
virtual void solveIter()
Execute one iteration of the solution algorithm.
virtual void solve()
Main control loop.
virtual void addPressureSource(fvScalarMatrix &matrix)
Source terms for the continuity equation.
objectiveManager objectiveManager_
Object to manage objective functions.
void allocateSensitivities()
Allocate the sensitivity derivatives.
virtual bool writeData(Ostream &os) const
Write the sensitivity derivatives.
autoPtr< adjointSensitivity > adjointSensitivity_
Sensitivity Derivatives engine.
const word & primalSolverName() const
Return the primal solver name.
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
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
const Type & value() const noexcept
Return const reference to value.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition fvMatrix.C:1415
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
tmp< volScalarField > A() const
Return the central coefficient.
Definition fvMatrix.C:1306
tmp< volScalarField > H1() const
Return H(1).
Definition fvMatrix.C:1377
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition fvMatrix.C:1004
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition fvMatrix.C:1260
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition fvMatrix.C:1325
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
List of finite volume options.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
static autoPtr< option > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvOption model.
Definition fvOption.C:78
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
const volScalarField & paInst() const
Return const reference to pressure.
Base class for incompressibleAdjoint solvers.
autoPtr< ATCModel > ATCModel_
Adjoint Transpose Convection options.
static autoPtr< incompressibleAdjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName)
Return a reference to the selected incompressible adjoint solver.
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
virtual bool includeDistance() const
Should the adjoint to the eikonal equation be solved.
virtual void updatePrimalBasedQuantities()
Update primal based quantities, e.g. the primal fields in adjoint turbulence models.
incompressibleVars & primalVars_
Primal variable set.
Class including all adjoint fields for incompressible flows.
autoPtr< variablesSet > vars_
Base variableSet pointer.
Definition solver.H:95
const word & managerType() const
Return the manager type.
Definition solverI.H:72
bool active_
Solve equations?
Definition solver.H:76
bool useSolverNameForFields() const
Use solver name as a suffix to the involved fields.
Definition solverI.H:36
const dictionary & dict() const
Return the solver dictionary.
Definition solverI.H:54
const word solverName_
Solver name.
Definition solver.H:71
const fvMesh & mesh() const
Return the solver mesh.
Definition solverI.H:24
const word & solverName() const
Return the solver name.
Definition solverI.H:30
fvMesh & mesh_
Reference to the mesh database.
Definition solver.H:56
virtual void postLoop()
Functions to be called after loop.
Definition solver.C:96
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const word dictName("faMeshDefinition")
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
auto & names
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition adjustPhi.C:30
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Type weightedAverage(const UList< scalar > &weights, const UList< Type > &fld)
The local weighted average of a field, using the mag() of the weights.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gMax(const FieldField< Field, Type > &f)
bool setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition findRefCell.C:27
fvScalarMatrix paEqn(fvm::d2dt2(pa) - sqr(c0) *fvc::laplacian(pa))
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
const bool momentumPredictor
dictionary dict
scalar sumLocalContErr
scalar globalContErr