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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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 "subCycleTime.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 namespace incompressible
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(adjointMeshMovementSolver, 0);
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  nLaplaceIters_ = dict_.getOrDefault<label>("iters", 1000);
50  tolerance_ = dict_.getOrDefault<scalar>("tolerance", 1e-6);
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 adjointMeshMovementSolver::adjointMeshMovementSolver
57 (
58  const fvMesh& mesh,
59  const dictionary& dict,
61  const labelHashSet& sensitivityPatchIDs,
62  const autoPtr<adjointEikonalSolver>& adjointEikonalSolverPtr
63 )
64 :
65  mesh_(mesh),
66  dict_(dict.subOrEmptyDict("adjointMeshMovementSolver")),
67  adjointSensitivity_(adjointSensitivity),
68  sensitivityPatchIDs_(sensitivityPatchIDs),
69  nLaplaceIters_(-1),
70  tolerance_(-1),
71  ma_
72  (
74  (
75  mesh,
76  "ma",
78  )
79  ),
80  source_
81  (
82  IOobject
83  (
84  "sourceAdjointMeshMovement",
85  mesh_.time().timeName(),
86  mesh_,
89  ),
90  mesh_,
92  ),
93  meshMovementSensPtr_(createZeroBoundaryPtr<vector>(mesh_)),
94  adjointEikonalSolverPtr_(adjointEikonalSolverPtr)
95 {
96  read();
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 {
104  dict_ = dict.subOrEmptyDict("adjointMeshMovementSolver");
105 
106  return true;
107 }
108 
109 
111 {
112  // Accumulate integrand from the current time step
114 
115  // Part of the source depending on the adjoint distance can be added only
116  // after solving the adjoint eikonal equation. Added in solve()
117 }
118 
119 
121 {
122  read();
123 
124  // Add source from the adjoint eikonal equation
126  {
127  source_ -=
128  fvc::div(adjointEikonalSolverPtr_().getFISensitivityTerm()().T());
129  }
130 
131  // Iterate the adjoint to the mesh movement equation
132  for (label iter = 0; iter < nLaplaceIters_; iter++)
133  {
134  Info<< "Adjoint Mesh Movement Iteration: " << iter << endl;
135 
136  fvVectorMatrix maEqn
137  (
139  + source_
140  );
141 
143 
144  //scalar residual = max(maEqn.solve().initialResidual());
145  scalar residual = mag(maEqn.solve().initialResidual());
146 
147  Info<< "Max ma " << gMax(mag(ma_)()) << endl;
148 
150 
151  // Check convergence
152  if (residual < tolerance_)
153  {
154  Info<< "\n***Reached adjoint mesh movement convergence limit, "
155  "iteration " << iter << "***\n\n";
156  break;
157  }
158  }
159  ma_.write();
160 }
161 
162 
164 {
165  source_ == dimensionedVector(source_.dimensions(), Zero);
167 }
168 
169 
171 {
172  Info<< "Calculating mesh movement sensitivities " << endl;
173 
174  boundaryVectorField& meshMovementSens = meshMovementSensPtr_();
175 
176  for (const label patchi : sensitivityPatchIDs_)
177  {
178  // No surface area included. Will be done by the actual sensitivity tool
179  meshMovementSens[patchi] = -ma_.boundaryField()[patchi].snGrad();
180  }
181 
182  return meshMovementSens;
183 }
184 
185 
187 {
188  return ma_;
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 } // End namespace incompressible
195 } // End namespace Foam
196 
197 // ************************************************************************* //
Foam::incompressible::adjointMeshMovementSolver::ma_
volVectorField ma_
Definition: adjointMeshMovementSolverIncompressible.H:80
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::fvMatrix::boundaryManipulate
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1348
Foam::incompressible::adjointMeshMovementSolver::adjointEikonalSolverPtr_
const autoPtr< adjointEikonalSolver > & adjointEikonalSolverPtr_
Definition: adjointMeshMovementSolverIncompressible.H:85
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::incompressible::adjointMeshMovementSolver::reset
void reset()
Reset source term.
Definition: adjointMeshMovementSolverIncompressible.C:163
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::incompressible::adjointMeshMovementSolver::accumulateIntegrand
void accumulateIntegrand(const scalar dt)
Accumulate source term.
Definition: adjointMeshMovementSolverIncompressible.C:110
Foam::incompressible::adjointMeshMovementSolver::read
void read()
Read options each time a new solution is found.
Definition: adjointMeshMovementSolverIncompressible.C:47
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::HashSet< label, Hash< label > >
Foam::incompressible::adjointMeshMovementSolver::solve
void solve()
Calculate the adjoint distance field.
Definition: adjointMeshMovementSolverIncompressible.C:120
Foam::incompressible::adjointMeshMovementSolver::source_
volVectorField source_
Definition: adjointMeshMovementSolverIncompressible.H:81
Foam::incompressible::adjointMeshMovementSolver::mesh_
const fvMesh & mesh_
Definition: adjointMeshMovementSolverIncompressible.H:74
Foam::incompressible::adjointMeshMovementSolver::meshMovementSensitivities
boundaryVectorField & meshMovementSensitivities()
Return the sensitivity term depending on da.
Definition: adjointMeshMovementSolverIncompressible.C:170
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::Time::printExecutionTime
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:618
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::incompressible::adjointSensitivity
Abstract base class for adjoint-based sensitivities in incompressible flows.
Definition: adjointSensitivityIncompressible.H:75
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::incompressible::adjointMeshMovementSolver::adjointSensitivity_
Foam::incompressible::adjointSensitivity & adjointSensitivity_
Definition: adjointMeshMovementSolverIncompressible.H:76
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::incompressible::adjointMeshMovementSolver::tolerance_
scalar tolerance_
Definition: adjointMeshMovementSolverIncompressible.H:79
Foam::incompressible::adjointMeshMovementSolver::nLaplaceIters_
label nLaplaceIters_
Definition: adjointMeshMovementSolverIncompressible.H:78
subCycleTime.H
Foam::incompressible::adjointMeshMovementSolver::readDict
virtual bool readDict(const dictionary &dict)
Read dict if changed.
Definition: adjointMeshMovementSolverIncompressible.C:102
Foam::incompressible::adjointMeshMovementSolver::ma
const volVectorField & ma()
Return the adjoint distance field.
Definition: adjointMeshMovementSolverIncompressible.C:186
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::variablesSet::autoCreateMeshMovementField
static tmp< volVectorField > autoCreateMeshMovementField(const fvMesh &mesh, const word &name, const dimensionSet &dims)
Auto create variable for mesh movement.
Definition: variablesSet.C:173
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressible::adjointMeshMovementSolver::meshMovementSensPtr_
autoPtr< boundaryVectorField > meshMovementSensPtr_
Wall face sens w.r.t.(x, y.z) //wall face sens w.r.t. (x,y.z)
Definition: adjointMeshMovementSolverIncompressible.H:84
Foam::incompressible::adjointMeshMovementSolver::sensitivityPatchIDs_
const labelHashSet & sensitivityPatchIDs_
Definition: adjointMeshMovementSolverIncompressible.H:77
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::incompressible::adjointMeshMovementSolver::dict_
dictionary dict_
Definition: adjointMeshMovementSolverIncompressible.H:75
Foam::GeometricField::Boundary
The boundary fields.
Definition: GeometricField.H:115
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::incompressible::adjointSensitivity::adjointMeshMovementSource
tmp< volVectorField > adjointMeshMovementSource()
Compute source term for adjoint mesh movement equation.
Definition: adjointSensitivityIncompressible.C:307
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
adjointMeshMovementSolverIncompressible.H