Loading...
Searching...
No Matches
adjointMeshMovementSolver.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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019 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
32#include "adjointSolver.H"
33#include "fvc.H"
34#include "fvm.H"
36#include "reverseLinear.H"
38
39// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40
41namespace Foam
42{
43
44defineTypeNameAndDebug(adjointMeshMovementSolver, 0);
45
46// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
47
49{
50 iters_ = dict_.getOrDefault<label>("iters", 1000);
51 tolerance_ = dict_.getOrDefault<scalar>("tolerance", 1.e-06);
52}
53
54
56{
58
59 // Add part related to the adjoint eikaonal equation, if necessary
60 const autoPtr<adjointEikonalSolver>& eikonalSolver =
62 if (eikonalSolver)
63 {
64 gradDxDbMult += eikonalSolver->getFISensitivityTerm();
65 }
66
67 source_ -=
69 (
70 mesh_.Sf()
71 & reverseLinear<tensor>(mesh_).interpolate(gradDxDbMult)
72 );
73
74 // Terms from objectives defined in (part of the) internal field
75 PtrList<objective>& functions =
77 getObjectiveFunctions();
78 for (objective& func : functions)
79 {
80 if (func.hasDivDxDbMult())
81 {
82 source_ -= func.weight()*fvc::grad(func.divDxDbMultiplier());
83 }
84 }
85
86 // Terms from fvOptions
87 source_.primitiveFieldRef() += adjointSensitivity_.optionsDxDbMult()();
88}
89
90
91// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
92
93adjointMeshMovementSolver::adjointMeshMovementSolver
94(
95 const fvMesh& mesh,
96 const dictionary& dict,
98)
99:
100 mesh_(mesh),
101 dict_(dict.subOrEmptyDict("adjointMeshMovementSolver")),
102 meshMovementSensPtr_(createZeroBoundaryPtr<vector>(mesh)),
103 adjointSensitivity_(adjointSensitivity),
104 ma_
105 (
106 variablesSet::autoCreateMeshMovementField
107 (
108 mesh_,
109 adjointSensitivity.getAdjointSolver().useSolverNameForFields()
110 ? ("ma" + adjointSensitivity.getAdjointSolver().solverName())
111 : "ma",
112 adjointSensitivity.getAdjointSolver().maDimensions()
113 )
114 ),
115 source_
116 (
118 (
119 "sourceadjointMeshMovement",
120 mesh_.time().timeName(),
121 mesh_,
122 IOobject::NO_READ,
123 IOobject::NO_WRITE
124 ),
125 mesh_,
127 (
128 adjointSensitivity.getAdjointSolver().maDimensions()/sqr(dimLength),
129 Zero
130 )
131 ),
132 iters_(0),
133 tolerance_(Zero)
135 read();
136}
137
138
139// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
140
142(
143 const dictionary& dict
144)
145{
146 dict_ = dict.subOrEmptyDict("adjointMeshMovementSolver");
147 read();
148
149 return true;
150}
151
152
154{
155 setSource();
156
157 // Iterate the adjoint to the mesh movement equation
158 for (label iter = 0; iter < iters_; iter++)
159 {
160 Info<< "adjoint Mesh Movement Iteration: " << iter << endl;
161
162 fvVectorMatrix maEqn
163 (
165 );
166
167 maEqn.boundaryManipulate(ma_.boundaryFieldRef());
168
169 scalar residual =
170 mag(Foam::solve(maEqn, mesh_.solverDict("ma")).initialResidual());
171
172 Info<< "Max ma " << gMax(mag(ma_)()) << endl;
173
175
176 // Check convergence
177 if (residual < tolerance_)
178 {
179 Info<< "\n***Reached adjoint mesh movement convergence limit, "
180 "iteration " << iter << "***\n\n";
181 break;
182 }
183 }
184 ma_.write();
185}
186
187
205 // No surface area included.
206 // Will be added during the assembly of the sensitivities
207 meshMovementSens[patchi] = -ma_.boundaryField()[patchi].snGrad();
208 }
209
210 return meshMovementSens;
211}
212
213
214// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215
216} // End namespace Foam
217
218// ************************************************************************* //
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Base class supporting Shape sensitivity derivatives.
virtual const labelHashSet & geometryVariationIntegrationPatches() const
Return set of patches on which to compute direct sensitivities.
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition TimeIO.C:608
Class solving the adjoint grid dispalcement PDEs. Assumes the primal grid displacement PDE is a Lapla...
volVectorField source_
Source term of the adjoint grid displacement PDEs.
autoPtr< boundaryVectorField > meshMovementSensPtr_
Part of sensitivity derivatives coming from the adjoint grid displacement PDE.
const fvMesh & mesh_
Reference to mesh.
ShapeSensitivitiesBase & adjointSensitivity_
dictionary dict_
Dictionary containing solution controls.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
boundaryVectorField & meshMovementSensitivities()
Return the sensitivity term depending on ma.
void read()
Read options each time a new solution is found.
volVectorField ma_
Adjoint grid displacement field.
void setSource()
Set the source term of the PDE.
virtual void solve()
Calculate the adjoint distance field.
Abstract base class for adjoint-based sensitivities.
autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver()
Return the adjoint eikonal solver.
const autoPtr< volTensorField > & gradDxDbMult() const
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
const objectiveManager & getObjectiveManager() const
Return a const reference to the objective manager.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition fvMatrix.C:1260
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const surfaceVectorField & Sf() const
Return cell face area vectors.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition objective.H:58
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Inversed weight central-differencing interpolation scheme class.
const dictionary & solverDict(const word &name) const
The solver controls dictionary for the given field. Same as solversDict().subDict(....
Definition solution.C:467
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
Base class for creating a set of variables.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Namespace for OpenFOAM.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
Vector< scalar > vector
Definition vector.H:57
volVectorField::Boundary boundaryVectorField
Type gMax(const FieldField< Field, Type > &f)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dictionary dict
Forwards and collection of common volume field types.