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 -------------------------------------------------------------------------------
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 
31 #include "wallFvPatch.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 namespace incompressible
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
50  wordList daTypes
51  (
52  mesh_.boundary().size(),
53  fixedValueFvPatchScalarField::typeName
54  );
55 
56  for (const label patchi : wallPatchIDs_)
57  {
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  (
85  IOobject
86  (
87  "ny",
88  mesh_.time().timeName(),
89  mesh_,
93  ),
94  mesh_,
96  patchDistMethod::patchTypes<vector>(mesh_, wallPatchIDs_)
97  );
98 
99  const fvPatchList& patches = mesh_.boundary();
100  volVectorField::Boundary& nybf = ny.boundaryFieldRef();
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 
117 adjointEikonalSolver::adjointEikonalSolver
118 (
119  const fvMesh& mesh,
120  const dictionary& dict,
121  const autoPtr<incompressible::RASModelVariables>& RASModelVars,
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  (
137  IOobject
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  (
156  IOobject
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 
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 {
240  distanceSensPtr_() = vector::zero;
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  (
274  IOobject
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  (
290  new volTensorField
291  (
292  IOobject
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 
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 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:218
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
autoPtr< boundaryVectorField > distanceSensPtr_
Wall face sens w.r.t. (x,y.z)
const surfaceVectorField & Sf() const
Return cell face area vectors.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
defineTypeNameAndDebug(adjointEikonalSolver, 0)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
const volScalarField & da()
Return the adjoint distance field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void solve()
Calculate the adjoint distance field.
wordList patchTypes(nPatches)
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:487
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
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:453
boundaryVectorField & distanceSensitivities()
Return the sensitivity term depending on da.
wordList patchTypes() const
Return the boundary condition types for da.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
Class including all adjoint fields for incompressible flows.
virtual tmp< volScalarField > distanceSensitivities()=0
Sensitivity terms resulting from the differentiation of the distance field. Misses dxdb...
word timeName
Definition: getTimeIndex.H:3
void read()
Read options each time a new solution is found.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
tmp< surfaceScalarField > computeYPhi()
Compute convecting velocity.
tmp< volVectorField > gradEikonal()
Return the gradient of the eikonal equation.
const autoPtr< incompressible::RASModelVariables > & RASModelVars_
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:46
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:113
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:212
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 current schemes dictionary, respects the "select" keyword.
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
void accumulateIntegrand(const scalar dt)
Accumulate source term.
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 time name of given scalar time formatted with the given precision.
Definition: Time.C:770
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)
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:620
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.
dimensionedScalar pow3(const dimensionedScalar &ds)
autoPtr< Foam::incompressibleAdjoint::adjointRASModel > & adjointTurbulence_
const volScalarField & d()
Return the distance field.
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:79
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:397
tmp< volTensorField > getFISensitivityTerm() const
Return the volume-based sensitivity term depending on da.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
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
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
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:133