volBSplinesBase.H
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 -------------------------------------------------------------------------------
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 Class
30  Foam::volBSplinesBase
31 
32 Description
33  Class constructing a number of volumetric B-Splines boxes,
34  read from dynamicMeshDict. Useful for various sensitivities and
35  optMeshMovement classes.
36 
37  Derives from MeshObject so that all instances know and update the same
38  control points and parametric coordinates are computed only once
39 
40 SourceFiles
41  volBSplinesBase.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef volBSplinesBase_H
46 #define volBSplinesBase_H
47 
48 #include "NURBS3DVolume.H"
49 #include "OFstream.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class volBSplinesBase Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 class volBSplinesBase
61 :
62  public MeshObject<fvMesh, UpdateableMeshObject, volBSplinesBase>
63 {
64 protected:
65 
66  // Protected data
67 
68  //- List with volumetric B-splines boxes.
69  // No overlapping is supported
71 
72  //- Active design variables numbering for all boxes
74 
75 
76 private:
77 
78  // Private Member Functions
79 
80  //- No copy construct
81  volBSplinesBase(const volBSplinesBase&) = delete;
82 
83  //- No copy assignment
84  void operator=(const volBSplinesBase&) = delete;
85 
86 
87 public:
88 
89  //- Runtime type information
90  TypeName("volBSplinesBase");
91 
92 
93  // Constructors
94 
95  //- Construct from components
96  volBSplinesBase(const fvMesh& mesh);
97 
98 
99  //- Destructor
100  virtual ~volBSplinesBase() = default;
101 
102 
103  // Member Functions
104 
105  //- Get const reference to the vol. B-splines boxes
106  const PtrList<NURBS3DVolume>& boxes() const;
107 
108  //- Get non-const reference to the vol. B-splines boxes
110 
111  //- Get const reference to a specific box
112  const NURBS3DVolume& box(const label boxI) const;
113 
114  //- Get non-const reference to a specific box
115  NURBS3DVolume& boxRef(const label boxI);
116 
117  //- Get reference to control points
118  const vectorField& getControlPoints(const label& iNURB) const;
119 
120  //- Get control points from all boxes
122 
123  //- Get cumulative number of control points from all boxes
124  label getTotalControlPointsNumber() const;
125 
126  //- Get number of boxes
127  label getNumberOfBoxes() const;
128 
129  //- Get start CP ID for each box
130  labelList getStartCpID() const;
131 
132  //- Get start CP ID for each box
133  labelList getStartVarID() const;
134 
135  //- Find box of certain control point
136  label findBoxID(const label cpI) const;
137 
138  //- From design variable ID, return boxID, cpID and direction
139  Vector<label> decomposeDV(const label dvI) const;
140 
141  //- Get active design variables
142  const labelList& getActiveDesignVariables() const;
143 
144  //- Get max boundary displacement for a given control-points
145  //- movement
147  (
148  const vectorField& controlPointsMovement,
149  const labelList& patchesToBeMoved
150  );
151 
152  //- Get the updated boundary points only
154  (
155  const vectorField& controlPointsMovement,
156  const labelList& patchesToBeMoved
157  );
158 
159  //- Bound control points movement
161  (
162  vectorField& controlPointsMovement
163  ) const;
164 
165  //- Move control points. No effect on mesh
166  void moveControlPoints(const vectorField& controlPointsMovement);
167 
168  //- Write control points to constant and optimisation folders
169  void writeControlPoints() const;
170 
171  //- Dummy function required by MeshObject.
172  // Since this class is going to initiate the mesh movement,
173  // there is nothing more to be done when the mesh points change
174  virtual bool movePoints();
175 
176  //- Dummy function required by MeshObject.
177  // Since this class is going to initiate the mesh movement,
178  // there is nothing more to be done when the mesh points change
179  virtual void updateMesh(const mapPolyMesh&);
180 };
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
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.
labelList getStartVarID() const
Get start CP ID for each box.
virtual void updateMesh(const mapPolyMesh &)
Dummy function required by MeshObject.
void writeControlPoints() const
Write control points to constant and optimisation folders.
labelList activeDesignVariables_
Active design variables numbering for all boxes.
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.
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
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:85
const fvMesh & mesh() const noexcept
Reference to the mesh.
Definition: MeshObject.H:157
labelList getStartCpID() const
Get start CP ID for each box.
const PtrList< NURBS3DVolume > & boxes() const
Get const reference to the vol. B-splines boxes.
TypeName("volBSplinesBase")
Runtime type information.
void boundControlPointMovement(vectorField &controlPointsMovement) const
Bound control points movement.
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:78
virtual ~volBSplinesBase()=default
Destructor.
virtual bool movePoints()
Dummy function required by MeshObject.
tmp< vectorField > computeBoundaryDisplacement(const vectorField &controlPointsMovement, const labelList &patchesToBeMoved)
Get the updated boundary points only.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
label findBoxID(const label cpI) const
Find box of certain control point.
PtrList< NURBS3DVolume > volume_
List with volumetric B-splines boxes.
Vector< label > decomposeDV(const label dvI) const
From design variable ID, return boxID, cpID and direction.
Namespace for OpenFOAM.
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.