volBSplinesBase.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 
30 #include "volBSplinesBase.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 defineTypeNameAndDebug(volBSplinesBase, 0);
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 volBSplinesBase::volBSplinesBase
44 (
45  const fvMesh& mesh
46 )
47 :
49  volume_(0),
50  activeDesignVariables_(0)
51 {
52  const dictionary NURBSdict
53  (
55  (
56  IOobject
57  (
58  "dynamicMeshDict",
59  mesh.time().constant(),
60  mesh,
63  false
64  )
65  ).subDict("volumetricBSplinesMotionSolverCoeffs")
66  );
67  // Read box names and allocate size
68  wordList controlBoxes(NURBSdict.toc());
69  volume_.setSize(controlBoxes.size());
70 
71  // Populate NURBS volumes
72  label iBox(0);
73  for (const word& boxName : controlBoxes)
74  {
75  if (NURBSdict.isDict(boxName))
76  {
77  volume_.set
78  (
79  iBox,
81  (
82  NURBSdict.subDict(boxName),
83  mesh,
84  true
85  )
86  );
87  volume_[iBox].write();
88  iBox++;
89  }
90  }
91  volume_.setSize(iBox);
92 
93  // Determine active design variables
94  activeDesignVariables_.setSize(3*getTotalControlPointsNumber(), -1);
95  label iActive(0);
96  const labelList startCpID(getStartCpID());
97  forAll(volume_, boxI)
98  {
99  const label start(3*startCpID[boxI]);
100  const boolList& isActiveVar = volume_[boxI].getActiveDesignVariables();
101  forAll(isActiveVar, varI)
102  {
103  if (isActiveVar[varI])
104  {
105  activeDesignVariables_[iActive++] = start + varI;
106  }
107  }
108  }
109  activeDesignVariables_.setSize(iActive);
110 }
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 {
117  return volume_;
118 }
119 
120 
122 {
123  return volume_;
124 }
125 
126 
127 const NURBS3DVolume& volBSplinesBase::box(const label boxI) const
128 {
129  return volume_[boxI];
130 }
131 
132 
134 {
135  return volume_[boxI];
136 }
137 
138 
139 const vectorField& volBSplinesBase::getControlPoints(const label& iNURB) const
140 {
141  return volume_[iNURB].getControlPoints();
142 }
143 
144 
146 {
147  vectorField totalCPs(0);
148  forAll(volume_, iNURB)
149  {
150  totalCPs.append(volume_[iNURB].getControlPoints());
151  }
152 
153  return totalCPs;
154 }
155 
156 
158 {
159  label nCPs(0);
160  forAll(volume_, iNURB)
161  {
162  nCPs += volume_[iNURB].getControlPoints().size();
163  }
164 
165  return nCPs;
166 }
167 
168 
170 {
171  return volume_.size();
172 }
173 
174 
176 {
177  // Allocate an extra entry to track in which box a CP might be
178  labelList startID(getNumberOfBoxes() + 1);
179  startID[0] = 0;
180  forAll(volume_, iNURB)
181  {
182  startID[iNURB+1] =
183  startID[iNURB] + volume_[iNURB].getControlPoints().size();
184  }
185 
186  return startID;
187 }
188 
189 
190 label volBSplinesBase::findBoxID(const label cpI) const
191 {
192  const labelList startCPID(getStartCpID());
193  for (label iBox = 0; iBox < startCPID.size() - 1 ; ++iBox)
194  {
195  if (cpI >= startCPID[iBox] || cpI < startCPID[iBox + 1])
196  {
197  return iBox;
198  }
199  }
200 
202  << "Invalid control point ID " << cpI << endl
203  << exit(FatalError);
204  return -1;
205 }
206 
207 
209 {
210  return activeDesignVariables_;
211 }
212 
213 
215 (
216  const vectorField& controlPointsMovement,
217  const labelList& patchesToBeMoved
218 )
219 {
220  scalar maxDisplacement(0);
221  label pastControlPoints(0);
222  forAll(volume_, iNURB)
223  {
224  const label nb(volume_[iNURB].getControlPoints().size());
225  vectorField localControlPointsMovement(nb, Zero);
226 
227  // Set localControlPointsMovement
228  forAll(localControlPointsMovement, iCPM)
229  {
230  localControlPointsMovement[iCPM] =
231  controlPointsMovement[pastControlPoints + iCPM];
232  }
233 
234  maxDisplacement = max
235  (
236  maxDisplacement,
237  volume_[iNURB].computeMaxBoundaryDisplacement
238  (
239  localControlPointsMovement,
240  patchesToBeMoved
241  )
242  );
243 
244  pastControlPoints += nb;
245  }
246 
247  return maxDisplacement;
248 }
249 
250 
252 (
253  vectorField& controlPointsMovement
254 )
255 {
256  label pastControlPoints(0);
257  forAll(volume_, iNURB)
258  {
259  const label nb(volume_[iNURB].getControlPoints().size());
260  vectorField localControlPointsMovement(nb, Zero);
261 
262  // Set localControlPointsMovement
263  forAll(localControlPointsMovement, iCPM)
264  {
265  localControlPointsMovement[iCPM] =
266  controlPointsMovement[pastControlPoints + iCPM];
267  }
268 
269  volume_[iNURB].boundControlPointMovement(localControlPointsMovement);
270 
271  // Transfer bounding back to controlPointMovement
272  forAll(localControlPointsMovement, iCPM)
273  {
274  controlPointsMovement[pastControlPoints + iCPM] =
275  localControlPointsMovement[iCPM];
276  }
277 
278  pastControlPoints += nb;
279  }
280 }
281 
282 
284 (
285  const vectorField& controlPointsMovement
286 )
287 {
288  label pastControlPoints(0);
289  forAll(volume_, iNURB)
290  {
291  const label nb(volume_[iNURB].getControlPoints().size());
292  vectorField localControlPointsMovement(nb, Zero);
293 
294  // Set localControlPointsMovement
295  forAll(localControlPointsMovement, iCPM)
296  {
297  localControlPointsMovement[iCPM] =
298  controlPointsMovement[pastControlPoints + iCPM];
299  }
300 
301  const vectorField newCps
302  (
303  volume_[iNURB].getControlPoints()
304  + localControlPointsMovement
305  );
306 
307  volume_[iNURB].setControlPoints(newCps);
308 
309  pastControlPoints += nb;
310  }
311 }
312 
313 
315 {
316  for (const NURBS3DVolume& box : volume_)
317  {
318  box.writeCps("cpsBsplines"+mesh_.time().timeName());
319  box.writeCpsInDict();
320  }
321 }
322 
323 
325 {
326  // Does nothing
327  return true;
328 }
329 
330 
332 {
333  // Does nothing
334 }
335 
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 
339 } // End namespace Foam
340 
341 // ************************************************************************* //
Foam::volBSplinesBase::activeDesignVariables_
labelList activeDesignVariables_
Active design variables numbering for all boxes.
Definition: volBSplinesBase.H:72
Foam::volBSplinesBase::movePoints
virtual bool movePoints()
Dummy function required by MeshObject.
Definition: volBSplinesBase.C:324
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::volBSplinesBase::boxRef
NURBS3DVolume & boxRef(const label boxI)
Get non-const reference to a specific box.
Definition: volBSplinesBase.C:133
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::volBSplinesBase::getAllControlPoints
vectorField getAllControlPoints() const
Get control points from all boxes.
Definition: volBSplinesBase.C:145
Foam::volBSplinesBase::findBoxID
label findBoxID(const label cpI) const
Find box of certain control point.
Definition: volBSplinesBase.C:190
Foam::NURBS3DVolume
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:72
Foam::volBSplinesBase::box
const NURBS3DVolume & box(const label boxI) const
Get const reference to a specific box.
Definition: volBSplinesBase.C:127
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::volBSplinesBase::getActiveDesignVariables
const labelList & getActiveDesignVariables() const
Get active design variables.
Definition: volBSplinesBase.C:208
Foam::volBSplinesBase::moveControlPoints
void moveControlPoints(const vectorField &controlPointsMovement)
Move control points. No effect on mesh.
Definition: volBSplinesBase.C:284
Foam::volBSplinesBase::boxes
const PtrList< NURBS3DVolume > & boxes() const
Get const reference to the vol. B-splines boxes.
Definition: volBSplinesBase.C:115
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Field< vector >
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::volBSplinesBase::getTotalControlPointsNumber
label getTotalControlPointsNumber() const
Get cumulative number of control points from all boxes.
Definition: volBSplinesBase.C:157
Foam::volBSplinesBase::getControlPoints
const vectorField & getControlPoints(const label &iNURB) const
Get reference to control points.
Definition: volBSplinesBase.C:139
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Foam::volBSplinesBase::boundControlPointMovement
void boundControlPointMovement(vectorField &controlPointsMovement)
Bound control points movement.
Definition: volBSplinesBase.C:252
Foam::volBSplinesBase::getNumberOfBoxes
label getNumberOfBoxes() const
Get number of boxes.
Definition: volBSplinesBase.C:169
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::volBSplinesBase::boxesRef
PtrList< NURBS3DVolume > & boxesRef()
Get non-const reference to the vol. B-splines boxes.
Definition: volBSplinesBase.C:121
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
volBSplinesBase.H
Foam::volBSplinesBase::volume_
PtrList< NURBS3DVolume > volume_
List with volumetric B-splines boxes.
Definition: volBSplinesBase.H:69
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::volBSplinesBase::computeMaxBoundaryDisplacement
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Definition: volBSplinesBase.C:215
Foam::volBSplinesBase::getStartCpID
labelList getStartCpID() const
Get start CP ID for each box.
Definition: volBSplinesBase.C:175
Foam::List< word >
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:186
Foam::volBSplinesBase::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Dummy function required by MeshObject.
Definition: volBSplinesBase.C:331
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::NURBS3DVolume::New
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
Definition: NURBS3DVolume.C:750
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::MeshObject< fvMesh, UpdateableMeshObject, volBSplinesBase >
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::volBSplinesBase::writeControlPoints
void writeControlPoints() const
Write control points to constant and optimisation folders.
Definition: volBSplinesBase.C:314