adjointEikonalSolver.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 -------------------------------------------------------------------------------
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 "adjointEikonalSolver.H"
31 #include "adjointSolver.H"
32 #include "fvc.H"
33 #include "fvm.H"
34 #include "surfaceInterpolation.H"
35 #include "volFieldsFwd.H"
36 #include "wallFvPatch.H"
37 #include "patchDistMethod.H"
38 #include "fvOptions.H"
40 #include "sensitivityTopO.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
55 {
56  wordList daTypes
57  (
58  mesh_.boundary().size(),
59  fixedValueFvPatchScalarField::typeName
60  );
61 
62  for (const label patchi : wallPatchIDs_)
63  {
65  }
66 
67  return daTypes;
68 }
69 
70 
72 {
73  nEikonalIters_ = dict_.getOrDefault<label>("iters", 1000);
74  tolerance_ = dict_.getOrDefault<scalar>("tolerance", 1e-6);
75  const scalar defaultEps =
76  mesh_.schemesDict().subDict("wallDist").
77  subOrEmptyDict("advectionDiffusionCoeffs").
78  getOrDefault<scalar>("epsilon", 0.1);
79  epsilon_ = dict_.getOrDefault<scalar>("epsilon", defaultEps);
80 }
81 
82 
84 {
85  // Primal distance field
87  const volScalarField& d = td();
88 
90  (
91  IOobject
92  (
93  "ny",
94  mesh_.time().timeName(),
95  mesh_,
99  ),
100  mesh_,
102  patchDistMethod::patchTypes<vector>(mesh_, wallPatchIDs_)
103  );
104 
105  const fvPatchList& patches = mesh_.boundary();
106  volVectorField::Boundary& nybf = ny.boundaryFieldRef();
107 
108  for (const label patchi : wallPatchIDs_)
109  {
110  nybf[patchi] == -patches[patchi].nf();
111  }
112 
113  ny = fvc::grad(d);
114 
116 
117  return tmp<surfaceScalarField>::New("yPhi", mesh_.Sf() & nf);
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
123 adjointEikonalSolver::adjointEikonalSolver
124 (
125  const fvMesh& mesh,
126  const dictionary& dict,
128  const labelHashSet& sensitivityPatchIDs
129 )
130 :
131  mesh_(mesh),
132  dict_(dict.subOrEmptyDict("adjointEikonalSolver")),
133  adjointSolver_(adjointSolver),
134  sensitivityPatchIDs_(sensitivityPatchIDs),
135  nEikonalIters_(-1),
136  tolerance_(-1),
137  epsilon_(Zero),
138  wallPatchIDs_(mesh_.boundaryMesh().findPatchIDs<wallPolyPatch>()),
139  da_
140  (
141  IOobject
142  (
143  word
144  (
145  adjointSolver.useSolverNameForFields() ?
146  "da" + adjointSolver.solverName() :
147  "da"
148  ),
149  mesh_.time().timeName(),
150  mesh_,
151  IOobject::READ_IF_PRESENT,
152  IOobject::AUTO_WRITE
153 // adjointVars.writeFields() ?
154 // IOobject::AUTO_WRITE : IOobject::NO_WRITE
155  ),
156  mesh_,
157  dimensionedScalar(adjointSolver.daDimensions() , Zero),
158  patchTypes()
159  ),
160  source_
161  (
162  IOobject
163  (
164  "sourceEikonal",
165  mesh_.time().timeName(),
166  mesh_,
167  IOobject::NO_READ,
168  IOobject::NO_WRITE
169  ),
170  mesh_,
172  ),
173  distanceSensPtr_(createZeroBoundaryPtr<vector>(mesh_))
174 {
175  read();
176 }
177 
178 
179 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180 
182 {
183  dict_ = dict.subOrEmptyDict("adjointEikonalSolver");
184  read();
185 
186  return true;
187 }
188 
189 
191 {
192  // Accumulate integrand from the current time step
194 }
195 
196 
198 {
199  read();
200 
201  // Primal distance field
203  const volScalarField& d = td();
204 
205  // Convecting flux
207  const surfaceScalarField& yPhi = tyPhi();
208 
209  volScalarField scaleDims
210  (
211  IOobject
212  (
213  "scaleDims",
214  mesh_.time().timeName(),
215  mesh_,
219  ),
220  mesh_,
221  dimensionedScalar("scaleDims", dimTime/dimLength, scalar(1))
222  );
223 
225 
226  // Iterate the adjoint to the eikonal equation
227  for (label iter = 0; iter < nEikonalIters_; ++iter)
228  {
229  read();
230 
231  Info<< "Adjoint Eikonal Iteration : " << iter << endl;
232 
233  fvScalarMatrix daEqn
234  (
235  2*fvm::div(-yPhi, da_)
238  + source_
239  ==
240  fvOptions(scaleDims, da_)
241  );
242 
243  daEqn.relax();
244  fvOptions.constrain(daEqn);
245  scalar residual = daEqn.solve().initialResidual();
246  fvOptions.correct(da_);
247 
248  Info<< "Max da " << gMax(mag(da_)()) << endl;
249 
251 
252  // Check convergence
253  if (residual < tolerance_)
254  {
255  Info<< "\n***Reached adjoint eikonal convergence limit, iteration "
256  << iter << "***\n\n";
257  break;
258  }
259  }
260  if (debug)
261  {
262  da_.write();
263  }
264 }
265 
266 
268 {
270  distanceSensPtr_() = vector::zero;
271 }
272 
273 
275 {
276  Info<< "Calculating distance sensitivities " << endl;
277 
278  boundaryVectorField& distanceSens = distanceSensPtr_();
279 
281  const volScalarField& d = td();
282 
283  for (const label patchi : sensitivityPatchIDs_)
284  {
285  vectorField nf(mesh_.boundary()[patchi].nf());
286 
287  // No surface area included. Will be done by the actual sensitivity tool
288  distanceSens[patchi] =
289  -2.*da_.boundaryField()[patchi]
290  *d.boundaryField()[patchi].snGrad()
291  *d.boundaryField()[patchi].snGrad()*nf;
292  }
293  return distanceSens;
294 }
295 
296 
297 tmp<volTensorField> adjointEikonalSolver::getFISensitivityTerm() const
298 {
299  Info<< "Calculating distance sensitivities " << endl;
300 
301  const tmp<volScalarField> td(adjointSolver_.yWall());
302  const volScalarField& d = td();
303  const volVectorField gradD(fvc::grad(d));
304 
305  auto gradDDa
306  (
308  (
309  IOobject
310  (
311  "gradDDa",
312  mesh_.time().timeName(),
313  mesh_,
316  ),
317  mesh_,
319  patchDistMethod::patchTypes<vector>(mesh_, wallPatchIDs_)
320  )
321  );
322  gradDDa.ref() = fvc::grad(d*da_);
323 
324  auto tdistanceSens
325  (
327  (
328  IOobject
329  (
330  "distanceSensFI",
331  mesh_.time().timeName(),
332  mesh_,
335  ),
336  mesh_,
339  )
340  );
341  volTensorField& distanceSens = tdistanceSens.ref();
342 
343  distanceSens =
344  - 2.*da_*gradD*gradD
345  - epsilon_*gradDDa*gradD
346  // grad(gradD) is symmetric theoretically but not numerically when
347  // computed with the Gauss divergence theorem. The following maintains
348  // exactly the same behaviour as the one before the re-factoring of
349  // sensitivities but avoid calling the tranpose operator.
351  distanceSens.correctBoundaryConditions();
352 
353  return tdistanceSens;
354 }
355 
356 
358 (
359  const word& designVarsName
360 ) const
361 {
363  const volScalarField& d = td();
364 
367 
370  (
371  tres.ref(), dSens, fvOptions, d.name(), designVarsName
372  );
373 
374  return tres;
375 }
376 
379 {
380  return da_;
381 }
382 
383 
385 {
387  const volScalarField& d = td();
388 
389  volVectorField gradD(fvc::grad(d));
390  return tmp<volVectorField>::New("gradEikonal", 2*gradD & fvc::grad(gradD));
391 }
392 
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 } // End namespace Foam
397 
398 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
tmp< scalarField > topologySensitivities(const word &designVarsName) const
Return sensitivity contribution to topology optimisation.
tmp< volVectorField > gradEikonal()
Return the gradient of the eikonal equation.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const surfaceVectorField & Sf() const
Return cell face area vectors.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
void accumulateIntegrand(const scalar dt)
Accumulate source term.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:88
boundaryVectorField & distanceSensitivities()
Return the sensitivity term depending on da.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
wordList patchTypes(nPatches)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const volScalarField & da()
Return the adjoint distance field.
Base class for adjoint solvers.
Definition: adjointSolver.H:53
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
void solve()
Calculate the adjoint distance field.
Finite-volume options.
Definition: fvOptions.H:51
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:40
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
word timeName
Definition: getTimeIndex.H:3
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
wordList patchTypes() const
Return the boundary condition types for da.
dynamicFvMesh & mesh
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
tmp< surfaceScalarField > computeYPhi()
Compute convecting velocity.
const labelHashSet & sensitivityPatchIDs_
A class for handling words, derived from Foam::string.
Definition: word.H:63
void reset()
Reset source term.
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:46
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
tmp< volTensorField > getFISensitivityTerm() const
Return the volume-based sensitivity term depending on da.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
const dictionary & schemesDict() const
The entire dictionary or the optional "select" sub-dictionary.
const volScalarField & d()
Return the distance field.
Calculation of adjoint based sensitivities for topology optimisation. This returns just the field par...
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
autoPtr< boundaryVectorField > distanceSensPtr_
Wall face sens w.r.t. (x,y.z)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
virtual tmp< volScalarField > adjointEikonalSource()
Return the source the adjoint eikonal equation.
defineTypeNameAndDebug(combustionModel, 0)
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:607
static void postProcessSens(scalarField &sens, scalarField &auxSens, fv::options &fvOptions, const word &fieldName, const word &designVariablesName)
Add part of the sensitivities coming from fvOptions.
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.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const polyBoundaryMesh & patches
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
void read()
Read options each time a new solution is found.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual tmp< volScalarField > yWall() const
Return the distance field, to be used in the solution of the adjoint eikonal PDE. ...
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Solver of the adjoint to the eikonal PDE.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
virtual bool readDict(const dictionary &dict)
Read dict if changed.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Do not request registration (bool: false)
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127