Loading...
Searching...
No Matches
incompressibleAdjointSolver.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-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 "wallFvPatch.H"
33#include "sensitivityTopO.H"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
45 (
48 adjointSolver
49 );
50}
51
52
53// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54
56(
57 const labelHashSet& sensitivityPatchIDs
58)
59{
60 if (hasBCdxdbMult_.bad())
61 {
62 const volVectorField& Ua = getAdjointVars().Ua();
63 hasBCdxdbMult_ = false;
64 for (const label patchI : sensitivityPatchIDs)
65 {
66 const fvPatchVectorField& Uab = Ua.boundaryField()[patchI];
68 {
71 (const_cast<fvPatchVectorField&>(Uab));
72 tmp<tensorField> dxdbMult = abc.dxdbMult();
73 if (dxdbMult)
74 {
75 hasBCdxdbMult_ = true;
76 break;
77 }
78 }
79 }
80 }
81 return hasBCdxdbMult_;
82}
83
84
85// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86
87Foam::incompressibleAdjointSolver::incompressibleAdjointSolver
88(
89 fvMesh& mesh,
90 const word& managerType,
91 const dictionary& dict,
92 const word& primalSolverName,
93 const word& solverName
94)
95:
96 adjointSolver(mesh, managerType, dict, primalSolverName, solverName),
97 primalVars_
98 (
99 mesh.lookupObjectRef<incompressiblePrimalSolver>(primalSolverName).
100 getIncoVars()
101 ),
102 ATCModel_(nullptr),
104{}
105
106
107// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
108
111(
112 fvMesh& mesh,
113 const word& managerType,
114 const dictionary& dict,
115 const word& primalSolverName,
116 const word& solverName
117)
118{
119 const word solverType(dict.get<word>("solver"));
120 auto* ctorPtr = dictionaryConstructorTable(solverType);
121
122 if (!ctorPtr)
123 {
125 (
126 dict,
127 "incompressibleAdjointSolver",
128 solverType,
129 *dictionaryConstructorTablePtr_
130 ) << exit(FatalIOError);
131 }
132
133 return
134 autoPtr<incompressibleAdjointSolver>
135 (
136 ctorPtr(mesh, managerType, dict, primalSolverName, solverName)
137 );
138}
139
140
141// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142
146 return primalVars_;
147}
148
149
150
153{
154 const incompressibleAdjointVars& adjointVars =
156 return adjointVars;
157}
158
159
162{
163 incompressibleAdjointVars& adjointVars =
165 return adjointVars;
166}
167
168
169
174}
175
178{
179 return getAdjointVars().adjointTurbulence()->includeDistance();
180}
181
186}
187
188
193
194
197{
198 return getAdjointVars().adjointTurbulence()->distanceSensitivities();
199}
200
201
206}
207
217 if (vars_)
218 {
219 getAdjointVars().adjointTurbulence()->setChangedPrimalSolution();
220 ATCModel_().updatePrimalBasedQuantities();
221 getAdjointVars().updatePrimalBasedQuantities();
222 }
223}
224
225
227(
228 volTensorField& gradDxDbMult,
229 const scalar dt
230)
231{
232 /*
233 addProfiling
234 (
235 incompressibleAdjointSolver,
236 "incompressibleAdjointSolver::accumulateGradDxDbMultiplier"
237 );
238 */
239 autoPtr<incompressibleAdjoint::adjointRASModel>& adjointRAS
240 (
241 getAdjointVars().adjointTurbulence()
242 );
243
244 const volScalarField& p = primalVars_.p();
245 const volVectorField& U = primalVars_.U();
246 const volScalarField& pa = getAdjointVars().pa();
247 const volVectorField& Ua = getAdjointVars().Ua();
248
249 // We only need to modify the boundaryField of gradU locally.
250 // If grad(U) is cached then
251 // a. The .ref() call fails since the tmp is initialised from a
252 // const ref
253 // b. we would be changing grad(U) for all other places in the code
254 // that need it
255 // So, always allocate new memory and avoid registering the new field
256 tmp<volTensorField> tgradU =
257 volTensorField::New("gradULocal", fvc::grad(U));
258 volTensorField& gradU = tgradU.ref();
259 volTensorField::Boundary& gradUbf = gradU.boundaryFieldRef();
260
261 // Explicitly correct the boundary gradient to get rid of
262 // the tangential component
263 forAll(mesh_.boundary(), patchI)
264 {
265 const fvPatch& patch = mesh_.boundary()[patchI];
266 if (isA<wallFvPatch>(patch))
267 {
268 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
269 gradUbf[patchI] = tnf*U.boundaryField()[patchI].snGrad();
270 }
271 }
272
273 tmp<volScalarField> tnuEff = adjointRAS->nuEff();
274 tmp<volSymmTensorField> stress = tnuEff()*twoSymm(gradU);
275 tmp<volTensorField> tgradUa = fvc::grad(Ua);
276
277 // Return field
278 auto tflowTerm
279 (
281 (
283 (
284 "flowTerm",
285 mesh_.time().timeName(),
286 mesh_,
289 ),
290 mesh_,
293 )
294 );
295 volTensorField& flowTerm = tflowTerm.ref();
296 // Term 3
297 flowTerm = - tnuEff*(gradU & twoSymm(tgradUa()));
298 // Term 4
299 flowTerm += fvc::grad(Ua & stress()) - (tgradUa & stress());
300
301 // Release memory
302 stress.clear();
303
304 // Compute dxdb multiplier
305 flowTerm +=
306 // Term 1, ATC
307 ATCModel_->getFISensitivityTerm()
308 // Term 2
309 - fvc::grad(p)*Ua;
310
311 // Term 5
312 flowTerm += pa*tgradU;
313
314 // Term 6, from the adjoint turbulence model
315 flowTerm += T(adjointRAS->FISensitivityTerm());
316
317 // Term 7, term from objective functions
318 PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
319
320 for (objective& objI : functions)
321 {
322 if (objI.hasGradDxDbMult())
323 {
324 flowTerm += objI.weight()*objI.gradDxDbMultiplier();
325 }
326 }
327
328 flowTerm.correctBoundaryConditions();
329
330 gradDxDbMult += flowTerm.T()*dt;
331 //profiling::writeNow();
332}
333
334
336(
337 autoPtr<scalarField>& divDxDbMult,
338 const scalar dt
339)
340{
341 PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
342 for (objective& func : functions)
343 {
344 if (func.hasDivDxDbMult())
345 {
346 divDxDbMult() +=
347 func.weight()*func.divDxDbMultiplier().primitiveField()*dt;
348 }
349 }
350}
351
352
354(
357 autoPtr<boundaryVectorField>& dxdbDirectMult,
358 autoPtr<pointBoundaryVectorField>& pointDxDbDirectMult,
359 const labelHashSet& sensitivityPatchIDs,
360 const scalar dt
361)
362{
363 PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
364 for (const label patchI : sensitivityPatchIDs)
365 {
366 const scalarField magSfDt(mesh_.boundary()[patchI].magSf()*dt);
367 for (objective& func : functions)
368 {
369 const scalar wei = func.weight();
370 if (func.hasdSdbMult())
371 {
372 dSfdbMult()[patchI] += wei*func.dSdbMultiplier(patchI)*dt;
373 }
374 if (func.hasdndbMult())
375 {
376 dnfdbMult()[patchI] += wei*func.dndbMultiplier(patchI)*magSfDt;
377 }
378 if (func.hasdxdbDirectMult())
379 {
380 dxdbDirectMult()[patchI] +=
381 wei*func.dxdbDirectMultiplier(patchI)*magSfDt;
382 }
383 }
384 }
385}
386
387
389(
391 const labelHashSet& sensitivityPatchIDs,
392 const scalar dt
393)
394{
395 if (!hasBCdxdbMult(sensitivityPatchIDs))
396 {
397 return;
398 }
399
400 // Grab references
401 const volVectorField& Ua = getAdjointVars().Ua();
402 const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointTurbulence =
403 getAdjointVars().adjointTurbulence();
404
405 // Fields needed to calculate adjoint sensitivities
407 turbVars = primalVars_.RASModelVariables();
408 const singlePhaseTransportModel& lamTransp = primalVars_.laminarTransport();
409 volScalarField nuEff(lamTransp.nu() + turbVars->nut());
410 tmp<volTensorField> tgradUa = fvc::grad(Ua);
411 const volTensorField::Boundary& gradUabf = tgradUa.cref().boundaryField();
412
413 // Avoid updating the event number to keep consistency with cases caching
414 // gradUa
415 auto& UaBoundary = getAdjointVars().Ua().boundaryFieldRef(false);
416 auto& nuEffBoundary = nuEff.boundaryField();
417
418 for (const label patchI : sensitivityPatchIDs)
419 {
420 fvPatchVectorField& Uab = UaBoundary[patchI];
422 {
423 tmp<tensorField> dxdbMult =
425 if (dxdbMult)
426 {
427 const fvPatch& patch = mesh_.boundary()[patchI];
428 tmp<vectorField> tnf = patch.nf();
429 const scalarField& magSf = patch.magSf();
430
431 tmp<vectorField> DvDbMult =
432 nuEffBoundary[patchI]
433 *(Uab.snGrad() + (gradUabf[patchI] & tnf))
434 // - (nf*pa.boundaryField()[patchI])
435 + adjointTurbulence().adjointMomentumBCSource()[patchI];
436 bcDxDbMult()[patchI] += (DvDbMult & dxdbMult())*magSf*dt;
437 }
438 }
439 }
440}
441
442
444(
445 vectorField& optionsDxDbMult,
446 const scalar dt
447)
448{
449 // Terms from fvOptions - missing contributions from turbulence models
450 const incompressibleAdjointVars& av = getAdjointVars();
451 vectorField temp(mesh_.nCells(), Zero);
453 (
454 temp, av.UaInst().name(), av.solverName()
455 );
456 optionsDxDbMult += temp*dt;
457 temp = Zero;
460 temp, av.paInst().name(), av.solverName()
461 );
462 optionsDxDbMult += temp*dt;
463}
464
465
467(
468 scalarField& betaMult,
469 const word& designVariablesName,
470 const scalar dt
471)
472{
473 const incompressibleAdjointVars& adjointVars = getAdjointVars();
474 const volVectorField& U = primalVars_.U();
475 const volVectorField& Ua = adjointVars.Ua();
476 const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointTurbulence =
477 adjointVars.adjointTurbulence();
479
480 // Term from the momentum equations
481 scalarField momSens((U.primitiveField() & Ua.primitiveField())*dt);
483 (betaMult, momSens, fvOptions, U.name(), designVariablesName);
484 if (debug > 2)
485 {
486 volScalarField IvSens
487 (
489 (
490 "IvSens" + solverName(),
491 mesh_.time().timeName(),
492 mesh_,
495 ),
496 mesh_,
498 );
499 IvSens.primitiveFieldRef() = momSens;
500 IvSens.write();
501 }
502
503 // Term from the turbulence model.
504 // Includes already the derivative of the interpolation function
505 betaMult +=
506 (adjointTurbulence->topologySensitivities(designVariablesName))*dt;
507
508 // Terms resulting directly from the objective function
509 PtrList<objective>& functions = objectiveManager_.getObjectiveFunctions();
510 for (objective& objI : functions)
511 {
512 const scalar weight(objI.weight());
513 if (objI.hasdJdb())
514 {
515 betaMult += weight*objI.dJdb()*dt;
516 }
517
518 if (objI.hasdJdbField())
519 {
520 SubField<scalar> betaSens(objI.dJdbField(), mesh_.nCells(), 0);
521 betaMult += weight*betaSens*dt;
522 }
523 }
524}
525
526
527// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
fv::options & fvOptions
static tmp< GeometricField< tensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< tensor >::calculatedType())
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field).
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< tensor, fvPatchField, volMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition SubField.H:58
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
bool bad() const noexcept
True if the Switch does not represent a valid enumeration.
Definition Switch.H:312
virtual tmp< Field< typename Foam::outerProduct< Foam::vector, Type >::type > > dxdbMult() const
Return contribution to sensitivity derivatives.
Base class for adjoint solvers.
objectiveManager objectiveManager_
Object to manage objective functions.
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
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
void postProcessSens(Field< Type > &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
Post process sensitivity field related to the fvOption.
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 volVectorField & Ua() const
Return const reference to velocity.
const volVectorField & UaInst() const
Return const reference to velocity.
const volScalarField & paInst() const
Return const reference to pressure.
Base class for incompressibleAdjoint solvers.
autoPtr< ATCModel > ATCModel_
Adjoint Transpose Convection options.
bool hasBCdxdbMult(const labelHashSet &sensitivityPatchIDs)
Compute, if necessary, and return the hasBCdxdbMult_ bool.
virtual tmp< volScalarField > yWall() const
Return the distance field, to be used in the solution of the adjoint eikonal PDE.
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.
const incompressibleVars & getPrimalVars() const
Access to the incompressible primal variables set.
virtual tmp< volScalarField > adjointEikonalSource()
Return the source the adjoint eikonal equation.
const autoPtr< ATCModel > & getATCModel() const
Access to the ATC model.
virtual dimensionSet daDimensions() const
Return the dimensions of the adjoint distance field.
virtual void accumulateGeometryVariationsMultipliers(autoPtr< boundaryVectorField > &dSfdbMult, autoPtr< boundaryVectorField > &dnfdbMult, autoPtr< boundaryVectorField > &dxdbDirectMult, autoPtr< pointBoundaryVectorField > &pointDxDirectDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Accumulate the multipliers of geometric quantities defined at the boundary, usually through an object...
virtual void accumulateDivDxDbMultiplier(autoPtr< scalarField > &divDxDbMult, const scalar dt)
Compute the multiplier for div(dxdb).
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.
virtual void accumulateGradDxDbMultiplier(volTensorField &gradDxDbMult, const scalar dt)
Compute the multiplier for grad(dxdb).
virtual void accumulateBCSensitivityIntegrand(autoPtr< boundaryVectorField > &bcDxDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Contributions from boundary functions that inlcude geometric aspects in them and change when the geom...
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
Switch hasBCdxdbMult_
Auxiliary bool to avoid a potentially expensive part of the sensitivity computation.
virtual void accumulateOptionsDxDbMultiplier(vectorField &optionsDxDbMult, const scalar dt)
Contributions from fvOptions that inlcude geometric aspects in them and change when the geometry is d...
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
incompressibleVars & primalVars_
Primal variable set.
Class including all adjoint fields for incompressible flows.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
Base class for primal incompressible solvers.
Base class for solution control classes.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition objective.H:58
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static void postProcessSens(scalarField &sens, scalarField &auxSens, fv::options &fvOptions, const word &fieldName, const word &designVariablesName)
Add part of the sensitivities coming from fvOptions.
A simple single-phase transport model based on viscosityModel.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
autoPtr< variablesSet > vars_
Base variableSet pointer.
Definition solver.H:95
const word & managerType() const
Return the manager type.
Definition solverI.H:72
const dictionary & dict() const
Return the solver dictionary.
Definition solverI.H:54
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
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Base class for creating a set of variables.
const word & solverName() const
Return solver name.
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
U
Definition pEqn.H:72
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
adjointBoundaryCondition< vector > adjointVectorBoundaryCondition
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< vector > fvPatchVectorField
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
Calculation of adjoint based sensitivities for topology optimisation. This returns just the field par...
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299