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,
150 )
151 :
152  SIBase
153  (
154  mesh,
155  dict,
156  primalVars,
157  adjointVars,
160  ),
161  volBSplinesBase_
162  (
164  ),
165 
166  flowSens_(0),
167  dSdbSens_(0),
168  dndbSens_(0),
169  dxdbDirectSens_(0),
170  bcSens_(0),
171 
172  derivativesFolder_("optimisation"/type() + "Derivatives")
173 {
174  // No boundary field pointers need to be allocated
175  const label nCPs(volBSplinesBase_.getTotalControlPointsNumber());
176  derivatives_ = scalarField(3*nCPs, Zero);
177  flowSens_ = vectorField(nCPs, Zero);
178  dSdbSens_ = vectorField(nCPs, Zero);
179  dndbSens_ = vectorField(nCPs, Zero);
180  dxdbDirectSens_ = vectorField(nCPs, Zero);
181  bcSens_ = vectorField(nCPs, Zero);
182 
183  // Create folder to store sensitivities
184  mkDir(derivativesFolder_);
185 }
186 
187 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189 
191 {
192  // Assemble the sensitivity map
193  // Solves for the post-processing equations and adds their contribution to
194  // the sensitivity map
196 
197  // Finalise sensitivities including dxFace/db
198  const boundaryVectorField& faceSens =
199  surfaceSensitivity_.getWallFaceSensVecBoundary();
200 
201  label passedCPs(0);
203  forAll(boxes, iNURB)
204  {
205  vectorField sens =
206  boxes[iNURB].computeControlPointSensitivities
207  (
208  faceSens,
209  sensitivityPatchIDs_.toc()
210  );
211  // Transfer to global list
212  forAll(sens, cpI)
213  {
214  flowSens_[passedCPs + cpI] = sens[cpI];
215  }
216  passedCPs += sens.size();
217  }
219 
220  // Contribution from objective function
221  // Note:
222  // includeObjectiveContribution has to be set to false (false by default)
223  // in surfaceSensitivity, in order to avoid computing this term twice.
224  // Optionally avoided altogether if includeObjectiveContribution is set to
225  // false for sensitivityVolBSplines
227 
229 
230  // Transform sensitivites to scalarField in order to cooperate with
231  // updateMethod
232  forAll(flowSens_, cpI)
233  {
234  derivatives_[3*cpI] =
235  flowSens_[cpI].x()
236  + dSdbSens_[cpI].x()
237  + dndbSens_[cpI].x()
238  + dxdbDirectSens_[cpI].x()
239  + bcSens_[cpI].x();
240  derivatives_[3*cpI + 1] =
241  flowSens_[cpI].y()
242  + dSdbSens_[cpI].y()
243  + dndbSens_[cpI].y()
244  + dxdbDirectSens_[cpI].y()
245  + bcSens_[cpI].y();
246  derivatives_[3*cpI + 2] =
247  flowSens_[cpI].z()
248  + dSdbSens_[cpI].z()
249  + dndbSens_[cpI].z()
250  + dxdbDirectSens_[cpI].z()
251  + bcSens_[cpI].z();
252  }
253 }
254 
255 
257 {
263 
265 }
266 
267 
268 void sensitivityVolBSplines::write(const word& baseName)
269 {
270  Info<< "Writing control point sensitivities to file" << endl;
271  // Write sensitivity map
272  SIBase::write(baseName);
273  // Write control point sensitivities
274  if (Pstream::master())
275  {
276  OFstream derivFile
277  (
279  baseName + adjointVars_.solverName() + mesh_.time().timeName()
280  );
281  unsigned int widthDV =
282  max(int(Foam::name(derivatives_.size()).size()), int(3));
283  unsigned int width = IOstream::defaultPrecision() + 7;
284  derivFile
285  << setw(widthDV) << "#cp" << " "
286  << setw(width) << "total::x"<< " "
287  << setw(width) << "total::y"<< " "
288  << setw(width) << "total::z"<< " "
289  << setw(width) << "flow::x" << " "
290  << setw(width) << "flow::y" << " "
291  << setw(width) << "flow::z" << " "
292  << setw(width) << "dSdb::x" << " "
293  << setw(width) << "dSdb::y" << " "
294  << setw(width) << "dSdb::z" << " "
295  << setw(width) << "dndb::x" << " "
296  << setw(width) << "dndb::y" << " "
297  << setw(width) << "dndb::z" << " "
298  << setw(width) << "dxdbDirect::x" << " "
299  << setw(width) << "dxdbDirect::y" << " "
300  << setw(width) << "dxdbDirect::z" << " "
301  << setw(width) << "dvdb::x" << " "
302  << setw(width) << "dvdb::y" << " "
303  << setw(width) << "dvdb::z" << endl;
304 
305  label passedCPs(0);
306  label lastActive(-1);
308  forAll(boxes, iNURB)
309  {
310  label nb = boxes[iNURB].getControlPoints().size();
311  const boolList& activeCPs = boxes[iNURB].getActiveCPs();
312  for (label iCP = 0; iCP < nb; iCP++)
313  {
314  if (activeCPs[iCP])
315  {
316  label globalCP = passedCPs + iCP;
317  if (globalCP!=lastActive + 1)
318  {
319  derivFile << "\n";
320  }
321  lastActive = globalCP;
322 
323  derivFile
324  << setw(widthDV) << globalCP << " "
325  << setw(width) << derivatives_[3*globalCP] << " "
326  << setw(width) << derivatives_[3*globalCP + 1] << " "
327  << setw(width) << derivatives_[3*globalCP + 2] << " "
328  << setw(width) << flowSens_[globalCP].x() << " "
329  << setw(width) << flowSens_[globalCP].y() << " "
330  << setw(width) << flowSens_[globalCP].z() << " "
331  << setw(width) << dSdbSens_[globalCP].x() << " "
332  << setw(width) << dSdbSens_[globalCP].y() << " "
333  << setw(width) << dSdbSens_[globalCP].z() << " "
334  << setw(width) << dndbSens_[globalCP].x() << " "
335  << setw(width) << dndbSens_[globalCP].y() << " "
336  << setw(width) << dndbSens_[globalCP].z() << " "
337  << setw(width) << dxdbDirectSens_[globalCP].x() << " "
338  << setw(width) << dxdbDirectSens_[globalCP].y() << " "
339  << setw(width) << dxdbDirectSens_[globalCP].z() << " "
340  << setw(width) << bcSens_[globalCP].x() << " "
341  << setw(width) << bcSens_[globalCP].y() << " "
342  << setw(width) << bcSens_[globalCP].z()
343  << endl;
344  }
345  }
346  passedCPs += nb;
347  }
348  }
349 }
350 
351 
352 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353 
354 } // End namespace incompressible
355 } // End namespace Foam
356 
357 // ************************************************************************* //
Foam::incompressible::adjointSensitivity::derivatives_
scalarField derivatives_
Definition: adjointSensitivityIncompressible.H:84
Foam::incompressible::sensitivityVolBSplines::assembleSensitivities
virtual void assembleSensitivities()
Assemble sensitivities.
Definition: sensitivityVolBSplinesIncompressible.C:190
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:617
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:62
Foam::incompressible::sensitivityVolBSplines::write
virtual void write(const word &baseName=word::null)
Write sensitivities to file.
Definition: sensitivityVolBSplinesIncompressible.C:268
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:256
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:785
Foam::incompressible::SIBase
Base class for Surface Integral-based sensitivity derivatives.
Definition: SIBaseIncompressible.H:57
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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::fv::optionAdjointList
Definition: fvOptionAdjointList.H:59
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:151
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 (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::incompressible::SIBase::write
virtual void write(const word &baseName=word::null)
Write sensitivity map.
Definition: SIBaseIncompressible.C:164
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:62
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:121
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:86
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
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:87
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:439
Foam::incompressible::SIBase::includeObjective_
bool includeObjective_
Whether to include direct sensitivities or not.
Definition: SIBaseIncompressible.H:72
Foam::incompressible::sensitivityVolBSplines::dSdbSens_
vectorField dSdbSens_
Term depending on delta(n dS)/delta b.
Definition: sensitivityVolBSplinesIncompressible.H:72
Foam::GeometricField::Boundary
Definition: GeometricField.H:114
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:333
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::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
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
fvOptionsAdjoint
fv::IOoptionListAdjoint fvOptionsAdjoint(mesh)
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