sensitivityVolBSplinesIncompressible.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(sensitivityVolBSplines, 0);
46 (
47  adjointSensitivity,
48  sensitivityVolBSplines,
49  dictionary
50 );
51 
52 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
53 
55 {
57  {
58  label passedCPs = 0;
60  forAll(boxes, iNURB)
61  {
62  label nb = boxes[iNURB].getControlPoints().size();
63  for (label cpI = 0; cpI < nb; cpI++)
64  {
65  vector dSdbSensCP(Zero);
66  vector dndbSensCP(Zero);
67  for (const label patchI : sensitivityPatchIDs_)
68  {
69  tensorField dSdb
70  (
71  boxes[iNURB].dndbBasedSensitivities(patchI, cpI)
72  );
73  dSdbSensCP += gSum(dSfdbMult_()[patchI] & dSdb);
74 
75  tensorField dndb
76  (
77  boxes[iNURB].dndbBasedSensitivities
78  (
79  patchI,
80  cpI,
81  false
82  )
83  );
84  dndbSensCP += gSum((dnfdbMult_()[patchI] & dndb));
85  }
86  dSdbSens_[passedCPs + cpI] = dSdbSensCP;
87  dndbSens_[passedCPs + cpI] = dndbSensCP;
88  }
89  passedCPs += nb;
90  }
93 
94  passedCPs = 0;
95  forAll(boxes, iNURB)
96  {
97  vectorField sensDxDbDirect =
98  boxes[iNURB].computeControlPointSensitivities
99  (
100  dxdbDirectMult_(),
101  sensitivityPatchIDs_.toc()
102  );
103 
104  // Transfer to global list
105  forAll(sensDxDbDirect, cpI)
106  {
107  dxdbDirectSens_[passedCPs + cpI] = sensDxDbDirect[cpI];
108  }
109  passedCPs += sensDxDbDirect.size();
110  }
112  }
113 }
114 
115 
117 {
118  label passedCPs = 0;
120  forAll(boxes, iNURB)
121  {
122  vectorField sensBcsDxDb =
123  boxes[iNURB].computeControlPointSensitivities
124  (
125  bcDxDbMult_(),
126  sensitivityPatchIDs_.toc()
127  );
128 
129  // Transfer to global list
130  forAll(sensBcsDxDb, cpI)
131  {
132  bcSens_[passedCPs + cpI] = sensBcsDxDb[cpI];
133  }
134  passedCPs += sensBcsDxDb.size();
135  }
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
142 sensitivityVolBSplines::sensitivityVolBSplines
143 (
144  const fvMesh& mesh,
145  const dictionary& dict,
146  incompressibleVars& primalVars,
147  incompressibleAdjointVars& adjointVars,
149 )
150 :
151  SIBase
152  (
153  mesh,
154  dict,
155  primalVars,
156  adjointVars,
158  ),
159  volBSplinesBase_
160  (
162  ),
163 
164  flowSens_(0),
165  dSdbSens_(0),
166  dndbSens_(0),
167  dxdbDirectSens_(0),
168  bcSens_(0),
169 
170  derivativesFolder_("optimisation"/type() + "Derivatives")
171 {
172  // No boundary field pointers need to be allocated
173  const label nCPs(volBSplinesBase_.getTotalControlPointsNumber());
174  derivatives_ = scalarField(3*nCPs, Zero);
175  flowSens_ = vectorField(nCPs, Zero);
176  dSdbSens_ = vectorField(nCPs, Zero);
177  dndbSens_ = vectorField(nCPs, Zero);
178  dxdbDirectSens_ = vectorField(nCPs, Zero);
179  bcSens_ = vectorField(nCPs, Zero);
180 
181  // Create folder to store sensitivities
182  mkDir(derivativesFolder_);
183 }
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
189 {
190  // Assemble the sensitivity map
191  // Solves for the post-processing equations and adds their contribution to
192  // the sensitivity map
194 
195  // Finalise sensitivities including dxFace/db
196  const boundaryVectorField& faceSens =
197  surfaceSensitivity_.getWallFaceSensVecBoundary();
198 
199  label passedCPs(0);
201  forAll(boxes, iNURB)
202  {
203  vectorField sens =
204  boxes[iNURB].computeControlPointSensitivities
205  (
206  faceSens,
207  sensitivityPatchIDs_.toc()
208  );
209  // Transfer to global list
210  forAll(sens, cpI)
211  {
212  flowSens_[passedCPs + cpI] = sens[cpI];
213  }
214  passedCPs += sens.size();
215  }
217 
218  // Contribution from objective function
219  // Note:
220  // includeObjectiveContribution has to be set to false (false by default)
221  // in surfaceSensitivity, in order to avoid computing this term twice.
222  // Optionally avoided altogether if includeObjectiveContribution is set to
223  // false for sensitivityVolBSplines
225 
227 
228  // Transform sensitivites to scalarField in order to cooperate with
229  // updateMethod
230  forAll(flowSens_, cpI)
231  {
232  derivatives_[3*cpI] =
233  flowSens_[cpI].x()
234  + dSdbSens_[cpI].x()
235  + dndbSens_[cpI].x()
236  + dxdbDirectSens_[cpI].x()
237  + bcSens_[cpI].x();
238  derivatives_[3*cpI + 1] =
239  flowSens_[cpI].y()
240  + dSdbSens_[cpI].y()
241  + dndbSens_[cpI].y()
242  + dxdbDirectSens_[cpI].y()
243  + bcSens_[cpI].y();
244  derivatives_[3*cpI + 2] =
245  flowSens_[cpI].z()
246  + dSdbSens_[cpI].z()
247  + dndbSens_[cpI].z()
248  + dxdbDirectSens_[cpI].z()
249  + bcSens_[cpI].z();
250  }
251 }
252 
253 
255 {
261 
263 }
264 
265 
266 void sensitivityVolBSplines::write(const word& baseName)
267 {
268  Info<< "Writing control point sensitivities to file" << endl;
269  // Write sensitivity map
270  SIBase::write(baseName);
271  // Write control point sensitivities
272  if (Pstream::master())
273  {
274  OFstream derivFile
275  (
277  baseName + adjointVars_.solverName() + mesh_.time().timeName()
278  );
279  unsigned int widthDV =
280  max(int(Foam::name(derivatives_.size()).size()), int(3));
281  unsigned int width = IOstream::defaultPrecision() + 7;
282  derivFile
283  << setw(widthDV) << "#cp" << " "
284  << setw(width) << "total::x"<< " "
285  << setw(width) << "total::y"<< " "
286  << setw(width) << "total::z"<< " "
287  << setw(width) << "flow::x" << " "
288  << setw(width) << "flow::y" << " "
289  << setw(width) << "flow::z" << " "
290  << setw(width) << "dSdb::x" << " "
291  << setw(width) << "dSdb::y" << " "
292  << setw(width) << "dSdb::z" << " "
293  << setw(width) << "dndb::x" << " "
294  << setw(width) << "dndb::y" << " "
295  << setw(width) << "dndb::z" << " "
296  << setw(width) << "dxdbDirect::x" << " "
297  << setw(width) << "dxdbDirect::y" << " "
298  << setw(width) << "dxdbDirect::z" << " "
299  << setw(width) << "dvdb::x" << " "
300  << setw(width) << "dvdb::y" << " "
301  << setw(width) << "dvdb::z" << endl;
302 
303  label passedCPs(0);
304  label lastActive(-1);
306  forAll(boxes, iNURB)
307  {
308  label nb = boxes[iNURB].getControlPoints().size();
309  const boolList& activeCPs = boxes[iNURB].getActiveCPs();
310  for (label iCP = 0; iCP < nb; iCP++)
311  {
312  if (activeCPs[iCP])
313  {
314  label globalCP = passedCPs + iCP;
315  if (globalCP!=lastActive + 1)
316  {
317  derivFile << "\n";
318  }
319  lastActive = globalCP;
320 
321  derivFile
322  << setw(widthDV) << globalCP << " "
323  << setw(width) << derivatives_[3*globalCP] << " "
324  << setw(width) << derivatives_[3*globalCP + 1] << " "
325  << setw(width) << derivatives_[3*globalCP + 2] << " "
326  << setw(width) << flowSens_[globalCP].x() << " "
327  << setw(width) << flowSens_[globalCP].y() << " "
328  << setw(width) << flowSens_[globalCP].z() << " "
329  << setw(width) << dSdbSens_[globalCP].x() << " "
330  << setw(width) << dSdbSens_[globalCP].y() << " "
331  << setw(width) << dSdbSens_[globalCP].z() << " "
332  << setw(width) << dndbSens_[globalCP].x() << " "
333  << setw(width) << dndbSens_[globalCP].y() << " "
334  << setw(width) << dndbSens_[globalCP].z() << " "
335  << setw(width) << dxdbDirectSens_[globalCP].x() << " "
336  << setw(width) << dxdbDirectSens_[globalCP].y() << " "
337  << setw(width) << dxdbDirectSens_[globalCP].z() << " "
338  << setw(width) << bcSens_[globalCP].x() << " "
339  << setw(width) << bcSens_[globalCP].y() << " "
340  << setw(width) << bcSens_[globalCP].z()
341  << endl;
342  }
343  }
344  passedCPs += nb;
345  }
346  }
347 }
348 
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 } // End namespace incompressible
353 } // End namespace Foam
354 
355 // ************************************************************************* //
Foam::incompressible::adjointSensitivity::derivatives_
scalarField derivatives_
Definition: adjointSensitivityIncompressible.H:83
Foam::incompressible::sensitivityVolBSplines::assembleSensitivities
virtual void assembleSensitivities()
Assemble sensitivities.
Definition: sensitivityVolBSplinesIncompressible.C:188
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::incompressible::sensitivityVolBSplines::computeBCContributions
void computeBCContributions()
Definition: sensitivityVolBSplinesIncompressible.C:116
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::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::incompressible::sensitivityVolBSplines::write
virtual void write(const word &baseName=word::null)
Write sensitivities to file.
Definition: sensitivityVolBSplinesIncompressible.C:266
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::volBSplinesBase
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
Definition: volBSplinesBase.H:59
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
sensitivityVolBSplinesIncompressible.H
Foam::incompressible::sensitivityVolBSplines::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: sensitivityVolBSplinesIncompressible.C:254
Foam::incompressible::sensitivityVolBSplines::volBSplinesBase_
volBSplinesBase & volBSplinesBase_
Reference to underlaying volumetric B-Splines morpher.
Definition: sensitivityVolBSplinesIncompressible.H:66
Foam::MeshObject< fvMesh, UpdateableMeshObject, volBSplinesBase >::New
static const volBSplinesBase & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
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::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::SIBase::surfaceSensitivity_
sensitivitySurface surfaceSensitivity_
Surface sensitivities.
Definition: SIBaseIncompressible.H:66
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::incompressible::sensitivityVolBSplines::computeObjectiveContributions
void computeObjectiveContributions()
Definition: sensitivityVolBSplinesIncompressible.C:54
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::incompressible::sensitivityVolBSplines::flowSens_
vectorField flowSens_
Flow related term.
Definition: sensitivityVolBSplinesIncompressible.H:69
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::incompressible::sensitivityVolBSplines::derivativesFolder_
fileName derivativesFolder_
Definition: sensitivityVolBSplinesIncompressible.H:84
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::incompressible::sensitivityVolBSplines::bcSens_
vectorField bcSens_
Term dependng on the differentiation of boundary conditions.
Definition: sensitivityVolBSplinesIncompressible.H:82
Foam::incompressible::SIBase::clearSensitivities
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
Definition: SIBaseIncompressible.C:148
Foam::incompressible::sensitivityVolBSplines::dxdbDirectSens_
vectorField dxdbDirectSens_
Definition: sensitivityVolBSplinesIncompressible.H:79
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
Foam::volBSplinesBase::boundControlPointMovement
void boundControlPointMovement(vectorField &controlPointsMovement)
Bound control points movement.
Definition: volBSplinesBase.C:252
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
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::incompressible::sensitivityVolBSplines::dndbSens_
vectorField dndbSens_
Term depending on delta (n)/delta b.
Definition: sensitivityVolBSplinesIncompressible.H:75
Foam::volBSplinesBase::boxesRef
PtrList< NURBS3DVolume > & boxesRef()
Get non-const reference to the vol. B-splines boxes.
Definition: volBSplinesBase.C:121
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::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::incompressible::sensitivityVolBSplines::dSdbSens_
vectorField dSdbSens_
Term depending on delta(n dS)/delta b.
Definition: sensitivityVolBSplinesIncompressible.H:72
Foam::GeometricField::Boundary
The boundary fields.
Definition: GeometricField.H:115
Foam::Vector< scalar >
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::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::incompressible::shapeSensitivities::dnfdbMult_
autoPtr< boundaryVectorField > dnfdbMult_
Definition: shapeSensitivitiesIncompressible.H:66
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::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54