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-------------------------------------------------------------------------------
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 "fvOptions.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39namespace incompressible
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
45
46
47// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48
50{
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
79(
80 const fvMesh& mesh,
81 const dictionary& dict,
83)
84:
86 gradDxDbMult_
87 (
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
111bool FIBase::readDict(const dictionary& dict)
112{
113 if (sensitivity::readDict(dict))
114 {
115 if (eikonalSolver_)
116 {
117 eikonalSolver_().readDict(dict);
118 }
119
120 return true;
121 }
122
123 return false;
124}
125
126
127void 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
144 fv::options::New(this->mesh_).postProcessSens
145 (
147 );
148
149 // Accumulate source for the adjoint to the eikonal equation
151 {
152 eikonalSolver_->accumulateIntegrand(dt);
153 }
154
155 // Accumulate direct sensitivities
157
158 // Accumulate sensitivities due to boundary conditions
160}
161
162
164{
168
170 {
171 eikonalSolver_->reset();
172 }
173
175}
176
177
178// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179
180} // End namespace incompressible
181} // End namespace Foam
182
183// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Base class for adjoint solvers.
Definition: adjointSolver.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for incompressibleAdjoint solvers.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Base class for Field Integral-based sensitivity derivatives.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
bool includeDistance_
Include distance variation in sens computation.
volTensorField gradDxDbMult_
grad(dx/db) multiplier
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
void read()
Read options and update solver pointers if necessary.
scalarField divDxDbMult_
div(dx/db) multiplier
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
vectorField optionsDxDbMult_
dx/db multiplier coming from fvOptions
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
virtual void accumulateBCSensitivityIntegrand(const scalar dt)
Accumulate sensitivities enamating from the boundary conditions.
virtual void accumulateDirectSensitivityIntegrand(const scalar dt)
Accumulate direct sensitivities.
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
const fvMesh & mesh_
Definition: sensitivity.H:69
dictionary dict_
Definition: sensitivity.H:70
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.C:57
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
void func(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict
A non-counting (dummy) refCount.
Definition: refCount.H:59