Loading...
Searching...
No Matches
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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
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
34namespace Foam
35{
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39defineTypeNameAndDebug(volBSplinesBase, 0);
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43volBSplinesBase::volBSplinesBase
44(
45 const fvMesh& mesh
46)
47:
48 MeshObject_type(mesh),
49 volume_(0),
51{
52 const dictionary NURBSdict
53 (
55 (
57 (
58 "dynamicMeshDict",
59 mesh.time().constant(),
60 mesh,
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 }
105}
106
107
108// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111{
112 return volume_;
113}
114
119}
120
122const NURBS3DVolume& volBSplinesBase::box(const label boxI) const
123{
124 return volume_[boxI];
125}
126
128NURBS3DVolume& volBSplinesBase::boxRef(const label boxI)
129{
130 return volume_[boxI];
131}
132
134const vectorField& volBSplinesBase::getControlPoints(const label& iNURB) const
135{
136 return volume_[iNURB].getControlPoints();
137}
138
139
141{
142 DynamicList<vector> totalCPs(0);
143 forAll(volume_, iNURB)
144 {
145 totalCPs.push_back(volume_[iNURB].getControlPoints());
146 }
147
148 return vectorField(std::move(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
186{
187 return 3*getStartCpID();
188}
189
190
191label volBSplinesBase::findBoxID(const label cpI) const
192{
193 const labelList startCPID(getStartCpID());
194 for (label iBox = 0; iBox < startCPID.size() - 1 ; ++iBox)
195 {
196 if (cpI >= startCPID[iBox] || cpI < startCPID[iBox + 1])
197 {
198 return iBox;
199 }
200 }
201
203 << "Invalid control point ID " << cpI << endl
204 << exit(FatalError);
205 return -1;
206}
207
208
209Vector<label> volBSplinesBase::decomposeDV(const label varID) const
210{
211 Vector<label> decomposed;
212 labelList startVarID = getStartVarID();
213 label boxID(-1);
214 for (label iBox = 0; iBox < startVarID.size() - 1 ; ++iBox)
215 {
216 if (varID >= startVarID[iBox] && varID < startVarID[iBox + 1])
217 {
218 boxID = iBox;
219 break;
220 }
221 }
222 const label localVarID = varID - startVarID[boxID];
223 decomposed.x() = boxID;
224 decomposed.y() = localVarID/3;
225 decomposed.z() = localVarID%3;
227 << "varID " << varID
228 << " belongs to box " << decomposed.x()
229 << " cpLocal " << decomposed.y()
230 << " dir " << decomposed.z()
231 << endl;
232 return decomposed;
233}
234
237{
239}
240
241
243(
244 const vectorField& controlPointsMovement,
245 const labelList& patchesToBeMoved
246)
247{
248 scalar maxDisplacement(0);
249 label pastControlPoints(0);
250 forAll(volume_, iNURB)
251 {
252 const label nb(volume_[iNURB].getControlPoints().size());
253 vectorField localControlPointsMovement(nb, Zero);
254
255 // Set localControlPointsMovement
256 forAll(localControlPointsMovement, iCPM)
257 {
258 localControlPointsMovement[iCPM] =
259 controlPointsMovement[pastControlPoints + iCPM];
260 }
261
262 maxDisplacement = max
263 (
264 maxDisplacement,
265 volume_[iNURB].computeMaxBoundaryDisplacement
266 (
267 localControlPointsMovement,
268 patchesToBeMoved
269 )
270 );
271
272 pastControlPoints += nb;
273 }
274
275 return maxDisplacement;
276}
277
278
280(
281 const vectorField& controlPointsMovement,
282 const labelList& patchesToBeMoved
283)
284{
285 auto tdisplacement(tmp<vectorField>::New(mesh_.nPoints(), Zero));
286 vectorField& displacement = tdisplacement.ref();
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 displacement +=
302 volume_[iNURB].computeNewBoundaryPoints
303 (
304 localControlPointsMovement,
305 patchesToBeMoved
306 )
307 - mesh_.points();
308
309 pastControlPoints += nb;
310 }
311
312 return tdisplacement;
313}
314
315
317(
318 vectorField& controlPointsMovement
319) const
320{
321 label pastControlPoints(0);
322 forAll(volume_, iNURB)
323 {
324 const label nb(volume_[iNURB].getControlPoints().size());
325 vectorField localControlPointsMovement(nb, Zero);
326
327 // Set localControlPointsMovement
328 forAll(localControlPointsMovement, iCPM)
329 {
330 localControlPointsMovement[iCPM] =
331 controlPointsMovement[pastControlPoints + iCPM];
332 }
333
334 volume_[iNURB].boundControlPointMovement(localControlPointsMovement);
335
336 // Transfer bounding back to controlPointMovement
337 forAll(localControlPointsMovement, iCPM)
338 {
339 controlPointsMovement[pastControlPoints + iCPM] =
340 localControlPointsMovement[iCPM];
342
343 pastControlPoints += nb;
344 }
345}
346
347
349(
350 const vectorField& controlPointsMovement
351)
352{
353 label pastControlPoints(0);
354 forAll(volume_, iNURB)
355 {
356 const label nb(volume_[iNURB].getControlPoints().size());
357 vectorField localControlPointsMovement(nb, Zero);
358
359 // Set localControlPointsMovement
360 forAll(localControlPointsMovement, iCPM)
361 {
362 localControlPointsMovement[iCPM] =
363 controlPointsMovement[pastControlPoints + iCPM];
364 }
365
366 const vectorField newCps
367 (
368 volume_[iNURB].getControlPoints()
369 + localControlPointsMovement
370 );
371
372 volume_[iNURB].setControlPoints(newCps);
373
374 pastControlPoints += nb;
375 }
376}
377
378
380{
381 for (const NURBS3DVolume& box : volume_)
382 {
383 box.writeCps("cpsBsplines" + mesh_.time().timeName());
384 }
385}
386
387
389{
390 // Does nothing
391 return true;
392}
393
394
396{
397 // Does nothing
398}
399
400
401// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
402
403} // End namespace Foam
404
405// ************************************************************************* //
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void push_back(const T &val)
Copy append an element to the end of this list.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void setSize(label n)
Alias for resize().
Definition List.H:536
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
static autoPtr< NURBS3DVolume > New(const dictionary &dict, const fvMesh &mesh, bool computeParamCoors=true)
Return a reference to the selected NURBS model.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
virtual bool movePoints()
Dummy function required by MeshObject.
void moveControlPoints(const vectorField &controlPointsMovement)
Move control points. No effect on mesh.
const vectorField & getControlPoints(const label &iNURB) const
Get reference to control points.
const PtrList< NURBS3DVolume > & boxes() const
Get const reference to the vol. B-splines boxes.
vectorField getAllControlPoints() const
Get control points from all boxes.
const labelList & getActiveDesignVariables() const
Get active design variables.
void writeControlPoints() const
Write control points to constant and optimisation folders.
scalar computeMaxBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Get max boundary displacement for a given control-points movement.
labelList getStartCpID() const
Get start CP ID for each box.
labelList getStartVarID() const
Get start CP ID for each box.
label findBoxID(const label cpI) const
Find box of certain control point.
NURBS3DVolume & boxRef(const label boxI)
Get non-const reference to a specific box.
PtrList< NURBS3DVolume > & boxesRef()
Get non-const reference to the vol. B-splines boxes.
labelList activeDesignVariables_
Active design variables numbering for all boxes.
const NURBS3DVolume & box(const label boxI) const
Get const reference to a specific box.
label getTotalControlPointsNumber() const
Get cumulative number of control points from all boxes.
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement.
PtrList< NURBS3DVolume > volume_
List with volumetric B-splines boxes.
virtual void updateMesh(const mapPolyMesh &)
Dummy function required by MeshObject.
tmp< vectorField > computeBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Get the updated boundary points only.
label getNumberOfBoxes() const
Get number of boxes.
Vector< label > decomposeDV(const label dvI) const
From design variable ID, return boxID, cpID and direction.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299