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 
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 
68  // Populate NURBS volumes
69  volume_.resize(NURBSdict.size());
70 
71  label iBox(0);
72 
73  for (const entry& dEntry : NURBSdict)
74  {
75  if (dEntry.isDict())
76  {
77  volume_.set
78  (
79  iBox,
80  NURBS3DVolume::New(dEntry.dict(), mesh, true)
81  );
82  volume_[iBox].writeParamCoordinates();
83  iBox++;
84  }
85  }
86  volume_.resize(iBox);
87 
88  // Determine active design variables
90  label iActive(0);
91  const labelList startCpID(getStartCpID());
92  forAll(volume_, boxI)
93  {
94  const label start(3*startCpID[boxI]);
95  const boolList& isActiveVar = volume_[boxI].getActiveDesignVariables();
96  forAll(isActiveVar, varI)
97  {
98  if (isActiveVar[varI])
99  {
100  activeDesignVariables_[iActive++] = start + varI;
101  }
102  }
103  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 {
112  return volume_;
113 }
114 
117 {
118  return volume_;
119 }
120 
122 const NURBS3DVolume& volBSplinesBase::box(const label boxI) const
123 {
124  return volume_[boxI];
125 }
126 
128 NURBS3DVolume& volBSplinesBase::boxRef(const label boxI)
129 {
130  return volume_[boxI];
131 }
132 
134 const vectorField& volBSplinesBase::getControlPoints(const label& iNURB) const
135 {
136  return volume_[iNURB].getControlPoints();
137 }
138 
139 
141 {
142  vectorField totalCPs(0);
143  forAll(volume_, iNURB)
144  {
145  totalCPs.append(volume_[iNURB].getControlPoints());
146  }
147 
148  return totalCPs;
149 }
150 
151 
153 {
154  label nCPs(0);
155  forAll(volume_, iNURB)
156  {
157  nCPs += volume_[iNURB].getControlPoints().size();
158  }
159 
160  return nCPs;
161 }
162 
165 {
166  return volume_.size();
167 }
168 
169 
171 {
172  // Allocate an extra entry to track in which box a CP might be
173  labelList startID(getNumberOfBoxes() + 1);
174  startID[0] = 0;
175  forAll(volume_, iNURB)
176  {
177  startID[iNURB+1] =
178  startID[iNURB] + volume_[iNURB].getControlPoints().size();
179  }
180 
181  return startID;
182 }
183 
184 
185 label volBSplinesBase::findBoxID(const label cpI) const
186 {
187  const labelList startCPID(getStartCpID());
188  for (label iBox = 0; iBox < startCPID.size() - 1 ; ++iBox)
189  {
190  if (cpI >= startCPID[iBox] || cpI < startCPID[iBox + 1])
191  {
192  return iBox;
193  }
194  }
195 
197  << "Invalid control point ID " << cpI << endl
198  << exit(FatalError);
199  return -1;
200 }
201 
202 
204 {
205  return activeDesignVariables_;
206 }
207 
208 
210 (
211  const vectorField& controlPointsMovement,
212  const labelList& patchesToBeMoved
213 )
214 {
215  scalar maxDisplacement(0);
216  label pastControlPoints(0);
217  forAll(volume_, iNURB)
218  {
219  const label nb(volume_[iNURB].getControlPoints().size());
220  vectorField localControlPointsMovement(nb, Zero);
221 
222  // Set localControlPointsMovement
223  forAll(localControlPointsMovement, iCPM)
224  {
225  localControlPointsMovement[iCPM] =
226  controlPointsMovement[pastControlPoints + iCPM];
227  }
228 
229  maxDisplacement = max
230  (
231  maxDisplacement,
232  volume_[iNURB].computeMaxBoundaryDisplacement
233  (
234  localControlPointsMovement,
235  patchesToBeMoved
236  )
237  );
238 
239  pastControlPoints += nb;
240  }
241 
242  return maxDisplacement;
243 }
244 
245 
247 (
248  vectorField& controlPointsMovement
249 ) const
250 {
251  label pastControlPoints(0);
252  forAll(volume_, iNURB)
253  {
254  const label nb(volume_[iNURB].getControlPoints().size());
255  vectorField localControlPointsMovement(nb, Zero);
256 
257  // Set localControlPointsMovement
258  forAll(localControlPointsMovement, iCPM)
259  {
260  localControlPointsMovement[iCPM] =
261  controlPointsMovement[pastControlPoints + iCPM];
262  }
263 
264  volume_[iNURB].boundControlPointMovement(localControlPointsMovement);
265 
266  // Transfer bounding back to controlPointMovement
267  forAll(localControlPointsMovement, iCPM)
268  {
269  controlPointsMovement[pastControlPoints + iCPM] =
270  localControlPointsMovement[iCPM];
271  }
273  pastControlPoints += nb;
274  }
275 }
276 
277 
279 (
280  const vectorField& controlPointsMovement
281 )
282 {
283  label pastControlPoints(0);
284  forAll(volume_, iNURB)
285  {
286  const label nb(volume_[iNURB].getControlPoints().size());
287  vectorField localControlPointsMovement(nb, Zero);
288 
289  // Set localControlPointsMovement
290  forAll(localControlPointsMovement, iCPM)
291  {
292  localControlPointsMovement[iCPM] =
293  controlPointsMovement[pastControlPoints + iCPM];
294  }
295 
296  const vectorField newCps
297  (
298  volume_[iNURB].getControlPoints()
299  + localControlPointsMovement
300  );
301 
302  volume_[iNURB].setControlPoints(newCps);
303 
304  pastControlPoints += nb;
305  }
306 }
307 
308 
310 {
311  for (const NURBS3DVolume& box : volume_)
312  {
313  box.writeCps("cpsBsplines" + mesh_.time().timeName());
314  }
315 }
316 
317 
319 {
320  // Does nothing
321  return true;
322 }
323 
324 
326 {
327  // Does nothing
328 }
329 
330 
331 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332 
333 } // End namespace Foam
334 
335 // ************************************************************************* //
const NURBS3DVolume & box(const label boxI) const
Get const reference to a specific box.
void moveControlPoints(const vectorField &controlPointsMovement)
Move control points. No effect on mesh.
virtual void updateMesh(const mapPolyMesh &)
Dummy function required by MeshObject.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
void writeControlPoints() const
Write control points to constant and optimisation folders.
labelList activeDesignVariables_
Active design variables numbering for all boxes.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const vectorField & getControlPoints(const label &iNURB) const
Get reference to control points.
label getNumberOfBoxes() const
Get number of boxes.
NURBS3DVolume & boxRef(const label boxI)
Get non-const reference to a specific box.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:453
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:69
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
vectorField getAllControlPoints() const
Get control points from all boxes.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:84
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dynamicFvMesh & mesh
Reading required, file watched for runTime modification.
labelList getStartCpID() const
Get start CP ID for each box.
const PtrList< NURBS3DVolume > & boxes() const
Get const reference to the vol. B-splines boxes.
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement.
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
defineTypeNameAndDebug(combustionModel, 0)
label getTotalControlPointsNumber() const
Get cumulative number of control points from all boxes.
PtrList< NURBS3DVolume > & boxesRef()
Get non-const reference to the vol. B-splines boxes.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual bool movePoints()
Dummy function required by MeshObject.
List< label > labelList
A List of labels.
Definition: List.H:62
label findBoxID(const label cpI) const
Find box of certain control point.
PtrList< NURBS3DVolume > volume_
List with volumetric B-splines boxes.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
List< bool > boolList
A List of bools.
Definition: List.H:60
Namespace for OpenFOAM.
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Get max boundary displacement for a given control-points movement.
const labelList & getActiveDesignVariables() const
Get active design variables.
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict. Useful for various sensitivities and optMeshMovement classes.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157