adjointEikonalSolverIncompressible.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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 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
31#include "wallFvPatch.H"
32#include "patchDistMethod.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39namespace incompressible
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
45
46// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47
49{
50 wordList daTypes
51 (
53 fixedValueFvPatchScalarField::typeName
54 );
55
56 for (const label patchi : wallPatchIDs_)
57 {
58 daTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
59 }
60
61 return daTypes;
62}
63
64
66{
67 nEikonalIters_ = dict_.getOrDefault<label>("iters", 1000);
68 tolerance_ = dict_.getOrDefault<scalar>("tolerance", 1e-6);
69 epsilon_ = dict_.getOrDefault<scalar>("epsilon", 0.1);
70 const scalar defaultEps =
71 mesh_.schemesDict().subDict("wallDist").
72 subOrEmptyDict("advectionDiffusionCoeffs").
73 getOrDefault<scalar>("epsilon", 0.1);
74 epsilon_ = dict_.getOrDefault<scalar>("epsilon", defaultEps);
75}
76
77
79{
80 // Primal distance field
81 const volScalarField& d = RASModelVars_().d();
82
84 (
86 (
87 "ny",
89 mesh_,
92 false
93 ),
94 mesh_,
96 patchDistMethod::patchTypes<vector>(mesh_, wallPatchIDs_)
97 );
98
101
102 for (const label patchi : wallPatchIDs_)
103 {
104 nybf[patchi] == -patches[patchi].nf();
105 }
106
107 ny = fvc::grad(d);
108
110
111 return tmp<surfaceScalarField>::New("yPhi", mesh_.Sf() & nf);
112}
113
114
115// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
116
118(
119 const fvMesh& mesh,
120 const dictionary& dict,
122 incompressibleAdjointVars& adjointVars,
123 const labelHashSet& sensitivityPatchIDs
124)
125:
126 mesh_(mesh),
127 dict_(dict.subOrEmptyDict("adjointEikonalSolver")),
128 RASModelVars_(RASModelVars),
129 adjointTurbulence_(adjointVars.adjointTurbulence()),
130 sensitivityPatchIDs_(sensitivityPatchIDs),
131 nEikonalIters_(-1),
132 tolerance_(-1),
133 epsilon_(Zero),
134 wallPatchIDs_(mesh_.boundaryMesh().findPatchIDs<wallPolyPatch>()),
135 da_
136 (
138 (
139 word
140 (
141 adjointVars.useSolverNameForFields() ?
142 "da" + adjointTurbulence_().adjointSolverName() :
143 "da"
144 ),
145 mesh_.time().timeName(),
146 mesh_,
147 IOobject::READ_IF_PRESENT,
148 IOobject::AUTO_WRITE
149 ),
150 mesh_,
152 patchTypes()
153 ),
154 source_
155 (
157 (
158 "sourceEikonal",
159 mesh_.time().timeName(),
160 mesh_,
161 IOobject::NO_READ,
162 IOobject::NO_WRITE
163 ),
164 mesh_,
166 ),
167 distanceSensPtr_(createZeroBoundaryPtr<vector>(mesh_))
168{
169 read();
170}
171
172
173// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
174
175bool adjointEikonalSolver::readDict(const dictionary& dict)
176{
177 dict_ = dict.subOrEmptyDict("adjointEikonalSolver");
178
179 return true;
180}
181
182
184{
185 // Accumulate integrand from the current time step
187}
188
189
191{
192 read();
193
194 // Primal distance field
195 const volScalarField& d = RASModelVars_().d();
196
197 // Convecting flux
199 const surfaceScalarField& yPhi = tyPhi();
200
201 // Iterate the adjoint to the eikonal equation
202 for (label iter = 0; iter < nEikonalIters_; ++iter)
203 {
204 read();
205
206 Info<< "Adjoint Eikonal Iteration : " << iter << endl;
207
208 fvScalarMatrix daEqn
209 (
210 2*fvm::div(-yPhi, da_)
213 + source_
214 );
215
216 daEqn.relax();
217 scalar residual = daEqn.solve().initialResidual();
218 Info<< "Max da " << gMax(mag(da_)()) << endl;
219
221
222 // Check convergence
223 if (residual < tolerance_)
224 {
225 Info<< "\n***Reached adjoint eikonal convergence limit, iteration "
226 << iter << "***\n\n";
227 break;
228 }
229 }
230 if (debug)
231 {
232 da_.write();
233 }
234}
235
236
238{
241}
242
243
245{
246 Info<< "Calculating distance sensitivities " << endl;
247
248 boundaryVectorField& distanceSens = distanceSensPtr_();
249
250 const volScalarField& d = RASModelVars_().d();
251 for (const label patchi : sensitivityPatchIDs_)
252 {
253 vectorField nf(mesh_.boundary()[patchi].nf());
254
255 // No surface area included. Will be done by the actual sensitivity tool
256 distanceSens[patchi] =
257 -2.*da_.boundaryField()[patchi]
258 *d.boundaryField()[patchi].snGrad()
259 *d.boundaryField()[patchi].snGrad()*nf;
260 }
261 return distanceSens;
262}
263
264
266{
267 Info<< "Calculating distance sensitivities " << endl;
268
269 const volScalarField& d = RASModelVars_().d();
270 const volVectorField gradD(fvc::grad(d));
271
272 volVectorField gradDDa
273 (
275 (
276 "gradDDa",
277 mesh_.time().timeName(),
278 mesh_,
281 ),
282 mesh_,
284 patchDistMethod::patchTypes<vector>(mesh_, wallPatchIDs_)
285 );
286 gradDDa = fvc::grad(d*da_);
287
288 tmp<volTensorField> tdistanceSens
289 (
291 (
293 (
294 "distanceSensFI",
295 mesh_.time().timeName(),
296 mesh_,
299 ),
300 mesh_,
302 )
303 );
304 volTensorField& distanceSens = tdistanceSens.ref();
305
306 distanceSens =
307 - 2.*da_*gradD*gradD
308 - epsilon_*gradD*gradDDa
309 + epsilon_*da_*d*fvc::grad(gradD);
310
311 return tdistanceSens;
312}
313
314
316{
317 return da_;
318}
319
320
322{
323 const volScalarField& d = RASModelVars_().d();
324 volVectorField gradD(fvc::grad(d));
325 return tmp<volVectorField>::New("gradEikonal", 2*gradD & fvc::grad(gradD));
326}
327
328
329// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330
331} // End namespace incompressible
332} // End namespace Foam
333
334// ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:618
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const surfaceVectorField & Sf() const
Return cell face area vectors.
Class including all adjoint fields for incompressible flows.
virtual tmp< volScalarField > distanceSensitivities()=0
tmp< volVectorField > gradEikonal()
Return the gradient of the eikonal equation.
boundaryVectorField & distanceSensitivities()
Return the sensitivity term depending on da.
const volScalarField & da()
Return the adjoint distance field.
autoPtr< boundaryVectorField > distanceSensPtr_
Wall face sens w.r.t. (x,y.z)
tmp< surfaceScalarField > computeYPhi()
Compute convecting velocity.
autoPtr< Foam::incompressibleAdjoint::adjointRASModel > & adjointTurbulence_
const autoPtr< incompressible::RASModelVariables > & RASModelVars_
const volScalarField & d()
Return the distance field.
tmp< volTensorField > getFISensitivityTerm() const
Return the volume-based sensitivity term depending on da.
void read()
Read options each time a new solution is found.
void accumulateIntegrand(const scalar dt)
Accumulate source term.
wordList patchTypes() const
Return the boundary condition types for da.
void solve()
Calculate the adjoint distance field.
virtual bool write(const bool valid=true) const
Write using setting from DB.
const dictionary & schemesDict() const
The current schemes dictionary, respects the "select" keyword.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:53
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
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
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:54
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
Type gMax(const FieldField< Field, Type > &f)
wordList patchTypes(nPatches)
dictionary dict
volScalarField & e
Definition: createFields.H:11
A non-counting (dummy) refCount.
Definition: refCount.H:59