FIBaseIncompressible.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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 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 "FIBaseIncompressible.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 namespace incompressible
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(FIBase, 0);
44 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
51  dict_.getOrDefault<bool>
52  (
53  "includeDistance",
54  adjointVars_.adjointTurbulence().ref().includeDistance()
55  );
56 
57  // Allocate distance solver if needed
58  if (includeDistance_ && eikonalSolver_.empty())
59  {
60  eikonalSolver_.reset
61  (
63  (
64  mesh_,
65  dict_,
68  sensitivityPatchIDs_
69  )
70  );
71  }
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
77 FIBase::FIBase
78 (
79  const fvMesh& mesh,
80  const dictionary& dict,
81  incompressibleVars& primalVars,
82  incompressibleAdjointVars& adjointVars,
85 )
86 :
88  (
89  mesh,
90  dict,
91  primalVars,
92  adjointVars,
95  ),
96  gradDxDbMult_
97  (
98  IOobject
99  (
100  "gradDxDbMult",
101  mesh_.time().timeName(),
102  mesh_,
105  ),
106  mesh_,
108  ),
109  divDxDbMult_(mesh_.nCells(), Zero),
110  optionsDxDbMult_(mesh_.nCells(), Zero),
111 
112  includeDistance_(false),
113  eikonalSolver_(nullptr)
114 {
115  read();
116 }
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
124  {
125  if (eikonalSolver_.valid())
126  {
127  eikonalSolver_().readDict(dict);
128  }
129 
130  return true;
131  }
132 
133  return false;
134 }
135 
136 
137 void FIBase::accumulateIntegrand(const scalar dt)
138 {
139  // Accumulate multiplier of grad(dxdb)
141 
142  // Accumulate multiplier of div(dxdb)
144  for (objective& func : functions)
145  {
146  divDxDbMult_ +=
147  func.weight()*func.divDxDbMultiplier().primitiveField()*dt;
148  }
149 
150  // Terms from fvOptions
151  for (fv::optionAdjoint& optionAdj : fvOptionsAdjoint_)
152  {
154  optionAdj.dxdbMult(adjointVars_)().primitiveField()*dt;
155  }
156 
157  // Accumulate source for the adjoint to the eikonal equation
158  if (includeDistance_)
159  {
160  eikonalSolver_->accumulateIntegrand(dt);
161  }
162 
163  // Accumulate direct sensitivities
165 
166  // Accumulate sensitivities due to boundary conditions
168 
169 }
170 
171 
173 {
175  divDxDbMult_ = Zero;
177 
178  if (includeDistance_)
179  {
180  eikonalSolver_->reset();
181  }
182 
184 }
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 } // End namespace incompressible
190 } // End namespace Foam
191 
192 // ************************************************************************* //
Foam::sensitivity::dict
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.C:57
Foam::incompressible::FIBase::optionsDxDbMult_
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
Definition: FIBaseIncompressible.H:72
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::objectiveManager
class for managing incompressible objective functions.
Definition: objectiveManager.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::incompressible::FIBase::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: FIBaseIncompressible.C:172
Foam::incompressible::adjointSensitivity::fvOptionsAdjoint_
fv::optionAdjointList & fvOptionsAdjoint_
Definition: adjointSensitivityIncompressible.H:88
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:107
Foam::incompressible::FIBase::gradDxDbMult_
volTensorField gradDxDbMult_
grad(dx/db) multiplier
Definition: FIBaseIncompressible.H:66
Foam::incompressible::shapeSensitivities::accumulateBCSensitivityIntegrand
virtual void accumulateBCSensitivityIntegrand(const scalar dt)
Accumulate sensitivities enamating from the boundary conditions.
Definition: shapeSensitivitiesIncompressible.C:66
Foam::incompressible::FIBase::read
void read()
Read options and update solver pointers if necessary.
Definition: FIBaseIncompressible.C:48
Foam::incompressible::FIBase::divDxDbMult_
scalarField divDxDbMult_
div(dx/db) multiplier
Definition: FIBaseIncompressible.H:69
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::fv::optionAdjointList
Definition: fvOptionAdjointList.H:59
Foam::incompressibleVars::RASModelVariables
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Definition: incompressibleVars.C:444
Foam::incompressible::FIBase::eikonalSolver_
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
Definition: FIBaseIncompressible.H:78
Foam::incompressible::shapeSensitivities
Definition: shapeSensitivitiesIncompressible.H:55
Foam::incompressible::adjointEikonalSolver
Solver of the adjoint to the eikonal PDE.
Definition: adjointEikonalSolverIncompressible.H:146
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::sensitivity::mesh_
const fvMesh & mesh_
Definition: sensitivity.H:69
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::sensitivity::readDict
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
Definition: sensitivity.C:63
Foam::fv::optionAdjoint
Similar to fv::option but with additional functionality to contribute to the sensitivity deriavtives.
Definition: fvOptionAdjoint.H:58
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::incompressible::FIBase::readDict
virtual bool readDict(const dictionary &dict)
Read dict if changed.
Definition: FIBaseIncompressible.C:121
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
Foam::func
void func(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
FIBaseIncompressible.H
Foam::incompressible::shapeSensitivities::accumulateDirectSensitivityIntegrand
virtual void accumulateDirectSensitivityIntegrand(const scalar dt)
Accumulate direct sensitivities.
Definition: shapeSensitivitiesIncompressible.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::incompressible::adjointSensitivity::adjointVars_
incompressibleAdjointVars & adjointVars_
Definition: adjointSensitivityIncompressible.H:86
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressible::adjointSensitivity::primalVars_
incompressibleVars & primalVars_
Definition: adjointSensitivityIncompressible.H:85
Foam::sensitivity::dict_
dictionary dict_
Definition: sensitivity.H:70
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::incompressible::FIBase::includeDistance_
bool includeDistance_
Include distance variation in sens computation.
Definition: FIBaseIncompressible.H:75
Foam::incompressibleAdjointVars::adjointTurbulence
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
Definition: incompressibleAdjointVars.C:70
Foam::incompressible::FIBase::accumulateIntegrand
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
Definition: FIBaseIncompressible.C:137
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::incompressible::shapeSensitivities::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: shapeSensitivitiesIncompressible.C:160
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::IOobject::NO_READ
Definition: IOobject.H:123
fvOptionsAdjoint
fv::IOoptionListAdjoint fvOptionsAdjoint(mesh)
Foam::objective
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:59
Foam::dimensionedTensor
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedTensor.H:52
Foam::incompressible::adjointSensitivity::computeGradDxDbMultiplier
tmp< volTensorField > computeGradDxDbMultiplier()
Definition: adjointSensitivityIncompressible.C:146
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::objectiveManager::getObjectiveFunctions
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Definition: objectiveManager.C:247
Foam::incompressible::adjointSensitivity::objectiveManager_
objectiveManager & objectiveManager_
Definition: adjointSensitivityIncompressible.H:87