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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 {
55  wordList daTypes
56  (
57  mesh_.boundary().size(),
58  fixedValueFvPatchScalarField::typeName
59  );
60 
61  for (const label patchi : wallPatchIDs_)
62  {
63  daTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
64  }
65 
66  return daTypes;
67 }
68 
69 
71 {
72  nEikonalIters_ = dict_.getOrDefault<label>("iters", 1000);
73  tolerance_ = dict_.getOrDefault<scalar>("tolerance", 1e-6);
74  const scalar defaultEps =
75  mesh_.schemesDict().subDict("wallDist").
76  subOrEmptyDict("advectionDiffusionCoeffs").
77  getOrDefault<scalar>("epsilon", 0.1);
78  epsilon_ = dict_.getOrDefault<scalar>("epsilon", defaultEps);
79 }
80 
81 
83 {
84  // Primal distance field
86 
88  (
89  IOobject
90  (
91  "ny",
92  mesh_.time().timeName(),
93  mesh_,
97  ),
98  mesh_,
100  patchDistMethod::patchTypes<vector>(mesh_, wallPatchIDs_)
101  );
102 
103  const fvPatchList& patches = mesh_.boundary();
104  volVectorField::Boundary& nybf = ny.boundaryFieldRef();
105 
106  for (const label patchi : wallPatchIDs_)
107  {
108  nybf[patchi] == -patches[patchi].nf();
109  }
110 
111  ny = fvc::grad(d);
112 
114 
115  return tmp<surfaceScalarField>::New("yPhi", mesh_.Sf() & nf);
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
120 
121 adjointEikonalSolver::adjointEikonalSolver
122 (
123  const fvMesh& mesh,
124  const dictionary& dict,
126  const labelHashSet& sensitivityPatchIDs
127 )
128 :
129  mesh_(mesh),
130  dict_(dict.subOrEmptyDict("adjointEikonalSolver")),
131  adjointSolver_(adjointSolver),
132  sensitivityPatchIDs_(sensitivityPatchIDs),
133  nEikonalIters_(-1),
134  tolerance_(-1),
135  epsilon_(Zero),
136  wallPatchIDs_(mesh_.boundaryMesh().findPatchIDs<wallPolyPatch>()),
137  da_
138  (
139  IOobject
140  (
141  word
142  (
143  adjointSolver.useSolverNameForFields() ?
144  "da" + adjointSolver.solverName() :
145  "da"
146  ),
147  mesh_.time().timeName(),
148  mesh_,
149  IOobject::READ_IF_PRESENT,
150  IOobject::AUTO_WRITE
151 // adjointVars.writeFields() ?
152 // IOobject::AUTO_WRITE : IOobject::NO_WRITE
153  ),
154  mesh_,
155  dimensionedScalar(adjointSolver.daDimensions() , Zero),
156  patchTypes()
157  ),
158  source_
159  (
160  IOobject
161  (
162  "sourceEikonal",
163  mesh_.time().timeName(),
164  mesh_,
165  IOobject::NO_READ,
166  IOobject::NO_WRITE
167  ),
168  mesh_,
170  ),
171  distanceSensPtr_(createZeroBoundaryPtr<vector>(mesh_))
172 {
173  read();
174 }
175 
176 
177 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178 
180 {
181  dict_ = dict.subOrEmptyDict("adjointEikonalSolver");
182  read();
183 
184  return true;
185 }
186 
187 
189 {
190  // Accumulate integrand from the current time step
192 }
193 
194 
196 {
197  read();
198 
199  // Primal distance field
201 
202  // Convecting flux
204  const surfaceScalarField& yPhi = tyPhi();
205 
206  volScalarField scaleDims
207  (
208  IOobject
209  (
210  "scaleDims",
211  mesh_.time().timeName(),
212  mesh_,
215  false
216  ),
217  mesh_,
218  dimensionedScalar("scaleDims", dimTime/dimLength, scalar(1))
219  );
220 
222 
223  // Iterate the adjoint to the eikonal equation
224  for (label iter = 0; iter < nEikonalIters_; ++iter)
225  {
226  read();
227 
228  Info<< "Adjoint Eikonal Iteration : " << iter << endl;
229 
230  fvScalarMatrix daEqn
231  (
232  2*fvm::div(-yPhi, da_)
235  + source_
236  ==
237  fvOptions(scaleDims, da_)
238  );
239 
240  daEqn.relax();
241  fvOptions.constrain(daEqn);
242  scalar residual = daEqn.solve().initialResidual();
243  fvOptions.correct(da_);
244 
245  Info<< "Max da " << gMax(mag(da_)()) << endl;
246 
248 
249  // Check convergence
250  if (residual < tolerance_)
251  {
252  Info<< "\n***Reached adjoint eikonal convergence limit, iteration "
253  << iter << "***\n\n";
254  break;
255  }
256  }
257  if (debug)
258  {
259  da_.write();
260  }
261 }
262 
263 
265 {
267  distanceSensPtr_() = vector::zero;
268 }
269 
270 
272 {
273  Info<< "Calculating distance sensitivities " << endl;
274 
275  boundaryVectorField& distanceSens = distanceSensPtr_();
276 
278  for (const label patchi : sensitivityPatchIDs_)
279  {
280  vectorField nf(mesh_.boundary()[patchi].nf());
281 
282  // No surface area included. Will be done by the actual sensitivity tool
283  distanceSens[patchi] =
284  -2.*da_.boundaryField()[patchi]
285  *d.boundaryField()[patchi].snGrad()
286  *d.boundaryField()[patchi].snGrad()*nf;
287  }
288  return distanceSens;
289 }
290 
291 
292 tmp<volTensorField> adjointEikonalSolver::getFISensitivityTerm() const
293 {
294  Info<< "Calculating distance sensitivities " << endl;
295 
297  const volVectorField gradD(fvc::grad(d));
298 
299  auto gradDDa
300  (
302  (
303  IOobject
304  (
305  "gradDDa",
306  mesh_.time().timeName(),
307  mesh_,
310  ),
311  mesh_,
313  patchDistMethod::patchTypes<vector>(mesh_, wallPatchIDs_)
314  )
315  );
316  gradDDa.ref() = fvc::grad(d*da_);
317 
318  auto tdistanceSens
319  (
321  (
322  IOobject
323  (
324  "distanceSensFI",
325  mesh_.time().timeName(),
326  mesh_,
329  ),
330  mesh_,
333  )
334  );
335  volTensorField& distanceSens = tdistanceSens.ref();
336 
337  distanceSens =
338  - 2.*da_*gradD*gradD
339  - epsilon_*gradDDa*gradD
340  // grad(gradD) is symmetric theoretically but not numerically when
341  // computed with the Gauss divergence theorem. The following maintains
342  // exactly the same behaviour as the one before the re-factoring of
343  // sensitivities but avoid calling the tranpose operator.
345  distanceSens.correctBoundaryConditions();
346 
347  return tdistanceSens;
348 }
349 
350 
352 (
353  const word& designVarsName
354 ) const
355 {
357 
360 
363  (
364  tres.ref(), dSens, fvOptions, d.name(), designVarsName
365  );
366 
367  return tres;
368 }
369 
372 {
373  return da_;
374 }
375 
376 
378 {
380  volVectorField gradD(fvc::grad(d));
381  return tmp<volVectorField>::New("gradEikonal", 2*gradD & fvc::grad(gradD));
382 }
383 
384 
385 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
386 
387 } // End namespace Foam
388 
389 // ************************************************************************* //
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
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
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:82
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:81
static const char *const typeName
Typename for Field.
Definition: Field.H:86
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:604
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:172
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