sensitivityBezierIncompressible.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 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 "IOmanip.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39namespace incompressible
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
46(
50);
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const fvMesh& mesh,
57 const dictionary& dict,
59)
60:
62 //Bezier_(mesh, dict), // AJH Read locally?
63 Bezier_(mesh, mesh.lookupObject<IOdictionary>("optimisationDict")),
64 sens_(Bezier_.nBezier(), Zero),
65 flowSens_(Bezier_.nBezier(), Zero),
66 dSdbSens_(Bezier_.nBezier(), Zero),
67 dndbSens_(Bezier_.nBezier(), Zero),
68 dxdbDirectSens_(Bezier_.nBezier(), Zero),
69 bcSens_(Bezier_.nBezier(), Zero),
70 derivativesFolder_("optimisation"/type() + "Derivatives")
71{
73 // Create folder to store sensitivities
75}
76
77
78// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79
81{
82 // Assemble the sensitivity map
83 // Solves for the post-processing equations and adds their contribution to
84 // the sensitivity map
86
87 forAll(sens_, iCP)
88 {
89 // Face-based summation. More robust since the result is independent of
90 // the number of processors (does not hold for a point-based summation)
91 for (const label patchI : sensitivityPatchIDs_)
92 {
93 // Interpolate parameterization info to faces
94 tmp<tensorField> tdxidXj = Bezier_.dxdbFace(patchI, iCP);
95 const tensorField& dxidXj = tdxidXj();
96
97 // Patch sensitivity map
98 const vectorField& patchSensMap =
99 surfaceSensitivity_.getWallFaceSensVecBoundary()[patchI];
100 flowSens_[iCP] += gSum(patchSensMap & dxidXj);
101
103 {
104 // Contribution from objective function
105 // term from delta( n dS ) / delta b and
106 // term from delta( n ) / delta b
107 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
108 tmp<tensorField> tdSdb
109 (
111 );
112 const tensorField& dSdb = tdSdb();
113 dSdbSens_[iCP] += gSum(dSfdbMult_()[patchI] & dSdb);
114
115 tmp<tensorField> tdndb
116 (
117 Bezier_.dndbBasedSensitivities(patchI, iCP, false)
118 );
119 const tensorField& dndb = tdndb();
120 dndbSens_[iCP] += gSum((dnfdbMult_()[patchI] & dndb));
121
122 // Contribution from objective function
123 // term from delta( x ) / delta b
124 // Only for objectives directly including
125 // x, like moments
126 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
127 dxdbDirectSens_[iCP] +=
128 gSum((dxdbDirectMult_()[patchI] & dxidXj));
129 }
130
131 // Sensitivities from boundary conditions
132 bcSens_[iCP] += gSum(bcDxDbMult_()[patchI] & dxidXj);
133 }
134 }
136
137 // Transform sensitivities to scalarField in order to cooperate with
138 // updateMethod
139 label nBezier = Bezier_.nBezier();
140 forAll(sens_, cpI)
141 {
142 derivatives_[cpI] = sens_[cpI].x();
143 derivatives_[cpI + nBezier] = sens_[cpI].y();
144 derivatives_[cpI + 2*nBezier] = sens_[cpI].z();
145 const boolList& confineXmovement = Bezier_.confineXmovement();
146 const boolList& confineYmovement = Bezier_.confineYmovement();
147 const boolList& confineZmovement = Bezier_.confineZmovement();
148 if (confineXmovement[cpI])
149 {
150 derivatives_[cpI] *= scalar(0);
151 flowSens_[cpI].x() = Zero;
152 dSdbSens_[cpI].x() = Zero;
153 dndbSens_[cpI].x() = Zero;
154 dxdbDirectSens_[cpI].x() = Zero;
155 bcSens_[cpI].x() = Zero;
156 }
157 if (confineYmovement[cpI])
158 {
159 derivatives_[cpI + nBezier] *= scalar(0);
160 flowSens_[cpI].y() = Zero;
161 dSdbSens_[cpI].y() = Zero;
162 dndbSens_[cpI].y() = Zero;
163 dxdbDirectSens_[cpI].y() = Zero;
164 bcSens_[cpI].y() = Zero;
165 }
166 if (confineZmovement[cpI])
167 {
168 derivatives_[cpI + 2*nBezier] *= scalar(0);
169 flowSens_[cpI].z() = Zero;
170 dSdbSens_[cpI].z() = Zero;
171 dndbSens_[cpI].z() = Zero;
172 dxdbDirectSens_[cpI].z() = Zero;
173 bcSens_[cpI].z() = Zero;
174 }
175 }
176}
177
178
180{
181 sens_ = Zero;
182 flowSens_ = Zero;
183 dSdbSens_ = Zero;
184 dndbSens_ = Zero;
186 bcSens_ = Zero;
187
189}
190
191
192void sensitivityBezier::write(const word& baseName)
193{
194 Info<< "Writing control point sensitivities to file" << endl;
195 // Write sensitivity map
196 SIBase::write(baseName);
197 // Write control point sensitivities
198 if (Pstream::master())
199 {
200 OFstream derivFile
201 (
203 baseName + adjointVars_.solverName() + mesh_.time().timeName()
204 );
205 unsigned int widthDV = max(int(name(sens_.size()).size()), int(3));
206 unsigned int width = IOstream::defaultPrecision() + 7;
207 derivFile
208 << setw(widthDV) << "#dv" << " "
209 << setw(width) << "total" << " "
210 << setw(width) << "flow" << " "
211 << setw(width) << "dSdb" << " "
212 << setw(width) << "dndb" << " "
213 << setw(width) << "dxdbDirect" << " "
214 << setw(width) << "dvdb" << endl;
215 label nDV = derivatives_.size();
216 label nBezier = Bezier_.nBezier();
217 const boolListList& confineMovement = Bezier_.confineMovement();
218 label lastActive(-1);
219 for (label iDV = 0; iDV < nDV; iDV++)
220 {
221 label iCP = iDV%nBezier;
222 label idir = iDV/nBezier;
223 if (!confineMovement[idir][iCP])
224 {
225 if (iDV!=lastActive + 1) derivFile << "\n";
226 lastActive = iDV;
227 derivFile
228 << setw(widthDV) << iDV << " "
229 << setw(width) << derivatives_[iDV] << " "
230 << setw(width) << flowSens_[iCP].component(idir) << " "
231 << setw(width) << dSdbSens_[iCP].component(idir) << " "
232 << setw(width) << dndbSens_[iCP].component(idir) << " "
233 << setw(width) << dxdbDirectSens_[iCP].component(idir) << " "
234 << setw(width) << bcSens_[iCP].component(idir)
235 << endl;
236 }
237 }
238 }
239}
240
241
242// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243
244} // End namespace incompressible
245} // End namespace Foam
246
247// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label nBezier() const
Number of Bezier control points.
Definition: Bezier.C:127
const boolList & confineZmovement() const
Confine z movement.
Definition: Bezier.C:151
const boolList & confineYmovement() const
Confine y movement.
Definition: Bezier.C:145
const boolListList & confineMovement() const
Info about confining movement in all directions.
Definition: Bezier.C:157
const boolList & confineXmovement() const
Confine x movement.
Definition: Bezier.C:139
tmp< tensorField > dxdbFace(const label patchI, const label cpI, bool useChainRule=true) const
dxdb tensor for a Bezier parameterized patch
Definition: Bezier.C:279
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool returnDimensionedNormalSens=true) const
Definition: Bezier.C:164
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:545
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
Output to file stream, using an OSstream.
Definition: OFstream.H:57
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
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
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Base class for incompressibleAdjoint solvers.
Base class for Surface Integral-based sensitivity derivatives.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
sensitivitySurface surfaceSensitivity_
Surface sensitivities.
bool includeObjective_
Whether to include direct sensitivities or not.
Abstract base class for adjoint-based sensitivities in incompressible flows.
Calculation of adjoint based sensitivities for Bezier control points.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
virtual void assembleSensitivities()
Assemble sensitivities.
virtual void assembleSensitivities()
Assemble sensitivities.
autoPtr< boundaryVectorField > dSfdbMult_
Fields related to direct sensitivities.
const fvMesh & mesh_
Definition: sensitivity.H:69
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:515
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333