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-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 
30 #include "FIBaseIncompressible.H"
32 #include "fvOptions.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 namespace incompressible
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
45 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
49 void FIBase::read()
50 {
52  dict_.getOrDefault<bool>
53  (
54  "includeDistance",
55  adjointVars_.adjointTurbulence().ref().includeDistance()
56  );
57 
58  // Allocate distance solver if needed
60  {
61  eikonalSolver_.reset
62  (
64  (
65  mesh_,
66  dict_,
69  sensitivityPatchIDs_
70  )
71  );
72  }
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 
78 FIBase::FIBase
79 (
80  const fvMesh& mesh,
81  const dictionary& dict,
83 )
84 :
85  shapeSensitivities(mesh, dict, adjointSolver),
86  gradDxDbMult_
87  (
88  IOobject
89  (
90  "gradDxDbMult",
91  mesh_.time().timeName(),
92  mesh_,
93  IOobject::NO_READ,
94  IOobject::NO_WRITE
95  ),
96  mesh_,
98  ),
99  divDxDbMult_(mesh_.nCells(), Zero),
100  optionsDxDbMult_(mesh_.nCells(), Zero),
101 
102  includeDistance_(false),
103  eikonalSolver_(nullptr)
104 {
105  read();
106 }
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
111 bool FIBase::readDict(const dictionary& dict)
112 {
114  {
115  if (eikonalSolver_)
116  {
117  eikonalSolver_().readDict(dict);
118  }
119 
120  return true;
121  }
122 
123  return false;
124 }
125 
126 
127 void FIBase::accumulateIntegrand(const scalar dt)
128 {
129  // Accumulate multiplier of grad(dxdb)
131 
132  // Accumulate multiplier of div(dxdb)
134  for (objective& func : functions)
135  {
136  if (func.hasDivDxDbMult())
137  {
138  divDxDbMult_ +=
139  func.weight()*func.divDxDbMultiplier().primitiveField()*dt;
140  }
141  }
142 
143  // Terms from fvOptions
145  (
147  );
148 
149  // Accumulate source for the adjoint to the eikonal equation
150  if (includeDistance_)
151  {
152  eikonalSolver_->accumulateIntegrand(dt);
153  }
154 
155  // Accumulate direct sensitivities
157 
158  // Accumulate sensitivities due to boundary conditions
160 }
161 
162 
164 {
166  divDxDbMult_ = Zero;
167  optionsDxDbMult_ = vector::zero;
168 
169  if (includeDistance_)
170  {
171  eikonalSolver_->reset();
172  }
173 
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 } // End namespace incompressible
181 } // End namespace Foam
182 
183 // ************************************************************************* //
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
const fvMesh & mesh_
Definition: sensitivity.H:65
dictionary dict
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
Definition: sensitivity.C:56
defineTypeNameAndDebug(adjointEikonalSolver, 0)
volTensorField gradDxDbMult_
grad(dx/db) multiplier
virtual void accumulateDirectSensitivityIntegrand(const scalar dt)
Accumulate direct sensitivities.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
void postProcessSens(Field< Type > &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
Post process sensitivity field related to the fvOption.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void read()
Read options and update solver pointers if necessary.
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
Base class for adjoint solvers.
Definition: adjointSolver.H:51
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Generic dimensioned Type class.
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
Base class for incompressibleAdjoint solvers.
scalarField divDxDbMult_
div(dx/db) multiplier
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.C:50
virtual bool readDict(const dictionary &dict)
Read dict if changed.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
const word & solverName() const
Return solver name.
Definition: variablesSet.C:77
bool includeDistance_
Include distance variation in sens computation.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:55
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
tmp< volTensorField > computeGradDxDbMultiplier()
Compute the volTensorField multiplying grad(dxdb) for the volume-based approach to compute shape sens...
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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 dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dictionary dict_
Definition: sensitivity.H:66
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...
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
Base class for Field Integral-based sensitivity derivatives.
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
virtual void accumulateBCSensitivityIntegrand(const scalar dt)
Accumulate sensitivities enamating from the boundary conditions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133