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