adjointMeshMovementSolverIncompressible.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-2022 PCOpt/NTUA
9  Copyright (C) 2013-2022 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 
33 #include "subCycleTime.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 namespace incompressible
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 {
51  nLaplaceIters_ = dict_.getOrDefault<label>("iters", 1000);
52  tolerance_ = dict_.getOrDefault<scalar>("tolerance", 1e-6);
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
58 adjointMeshMovementSolver::adjointMeshMovementSolver
59 (
60  const fvMesh& mesh,
61  const dictionary& dict,
63  const labelHashSet& sensitivityPatchIDs,
64  const autoPtr<adjointEikonalSolver>& adjointEikonalSolverPtr
65 )
66 :
67  mesh_(mesh),
68  dict_(dict.subOrEmptyDict("adjointMeshMovementSolver")),
69  adjointSensitivity_(adjointSensitivity),
70  sensitivityPatchIDs_(sensitivityPatchIDs),
71  nLaplaceIters_(-1),
72  tolerance_(-1),
73  ma_
74  (
75  IOobject
76  (
77  word
78  (
79  adjointSensitivity.adjointVars().useSolverNameForFields()
80  ? "ma" + adjointSensitivity.adjointSolver().solverName()
81  : "ma"
82  ),
83  mesh.time().timeName(),
84  mesh,
85  IOobject::READ_IF_PRESENT,
86  IOobject::AUTO_WRITE
87  ),
88  mesh,
90  fixedValueFvPatchVectorField::typeName
91  ),
92  source_
93  (
94  IOobject
95  (
96  "sourceAdjointMeshMovement",
97  mesh_.time().timeName(),
98  mesh_,
99  IOobject::NO_READ,
100  IOobject::NO_WRITE
101  ),
102  mesh_,
104  ),
105  meshMovementSensPtr_(createZeroBoundaryPtr<vector>(mesh_)),
106  adjointEikonalSolverPtr_(adjointEikonalSolverPtr)
107 {
108  read();
109 }
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  dict_ = dict.subOrEmptyDict("adjointMeshMovementSolver");
117 
118  return true;
119 }
120 
121 
123 {
124  // Accumulate integrand from the current time step
126 
127  // Part of the source depending on the adjoint distance can be added only
128  // after solving the adjoint eikonal equation. Added in solve()
129 }
130 
131 
133 {
134  read();
135 
136  // Add source from the adjoint eikonal equation
138  {
139  source_ -=
140  fvc::div(adjointEikonalSolverPtr_().getFISensitivityTerm()().T());
141  }
142 
143  // Iterate the adjoint to the mesh movement equation
144  for (label iter = 0; iter < nLaplaceIters_; iter++)
145  {
146  Info<< "Adjoint Mesh Movement Iteration: " << iter << endl;
147 
148  fvVectorMatrix maEqn
149  (
151  + source_
152  );
153 
154  maEqn.boundaryManipulate(ma_.boundaryFieldRef());
155 
156  //scalar residual = max(maEqn.solve().initialResidual());
157  scalar residual =
158  mag(Foam::solve(maEqn, mesh_.solverDict("ma")).initialResidual());
159 
160  Info<< "Max ma " << gMax(mag(ma_)()) << endl;
161 
163 
164  // Check convergence
165  if (residual < tolerance_)
166  {
167  Info<< "\n***Reached adjoint mesh movement convergence limit, "
168  "iteration " << iter << "***\n\n";
169  break;
170  }
171  }
172  ma_.write();
173 }
174 
175 
177 {
179  meshMovementSensPtr_() = vector::zero;
180 }
181 
182 
184 {
185  Info<< "Calculating mesh movement sensitivities " << endl;
186 
187  boundaryVectorField& meshMovementSens = meshMovementSensPtr_();
188 
189  for (const label patchi : sensitivityPatchIDs_)
190  {
191  // No surface area included. Will be done by the actual sensitivity tool
192  meshMovementSens[patchi] = -ma_.boundaryField()[patchi].snGrad();
193  }
194 
195  return meshMovementSens;
196 }
197 
198 
200 {
201  return ma_;
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace incompressible
208 } // End namespace Foam
209 
210 // ************************************************************************* //
Solver of the adjoint to the Laplace grid displacement equation.
dictionary dict
defineTypeNameAndDebug(adjointEikonalSolver, 0)
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
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
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
void accumulateIntegrand(const scalar dt)
Accumulate source term.
autoPtr< boundaryVectorField > meshMovementSensPtr_
Wall face sens w.r.t.(x, y.z) //wall face sens w.r.t. (x,y.z)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void read()
Read options each time a new solution is found.
Base class for adjoint solvers.
Definition: adjointSolver.H:51
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:363
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
word timeName
Definition: getTimeIndex.H:3
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)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Abstract base class for adjoint-based sensitivities in incompressible flows.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
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
virtual bool readDict(const dictionary &dict)
Read dict if changed.
tmp< volVectorField > adjointMeshMovementSource()
Compute source term for adjoint mesh movement equation.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Type gMax(const FieldField< Field, Type > &f)
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:620
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
boundaryVectorField & meshMovementSensitivities()
Return the sensitivity term depending on da.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
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...
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.
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133
const volVectorField & ma()
Return the adjoint distance field.