optMeshMovementVolumetricBSplines.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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 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
28Class
29 Foam::optMeshMovementVolumetricBSplines
30
31Description
32 Converts NURBS volume control points update to actual mesh movement.
33 Internal points are also moved based on the movement of the control points
34
35SourceFiles
36 optMeshMovementVolumetricBSplines.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef optMeshMovementVolumetricBSplines_H
41#define optMeshMovementVolumetricBSplines_H
42
43#include "optMeshMovement.H"
44#include "volBSplinesBase.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class optMeshMovementVolumetricBSplines Declaration
53\*---------------------------------------------------------------------------*/
56:
57 public optMeshMovement
58{
59protected:
60
61 // Protected data
62
63 //- Reference to underlaying volumetric B-Splines morpher
65
66 //- Backup of initial control points. Useful for line-search
68
69
70 // Protected Member Functions
71
73
74
75private:
76
77 // Private Member Functions
78
79 //- No copy construct
81 (
83 ) = delete;
84
85 //- No copy assignment
86 void operator=(const optMeshMovementVolumetricBSplines&) = delete;
87
88
89public:
90
91 //- Runtime type information
92 TypeName("volumetricBSplines");
93
94
95 // Constructors
96
97 //- Construct from components
99 (
100 fvMesh& mesh,
101 const dictionary& dict,
102 const labelList& patchIDs
103 );
104
105
106 //- Destructor
107 virtual ~optMeshMovementVolumetricBSplines() = default;
108
109
110 // Member Functions
111
112 //- Calculates surface mesh movement
113 void moveMesh();
114
115 //- Store design variables and mesh, to act as the starting point of
116 //- line search
117 virtual void storeDesignVariables();
118
119 //- Reset to starting point of line search
120 virtual void resetDesignVariables();
121
122 //- Compute eta value based on max displacement
123 virtual scalar computeEta(const scalarField& correction);
124
125 //- Return active design variables
126 virtual labelList getActiveDesignVariables() const;
127};
128
129
130// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131
132} // End namespace Foam
133
134// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135
136#endif
137
138// ************************************************************************* //
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Converts NURBS volume control points update to actual mesh movement. Internal points are also moved b...
virtual void resetDesignVariables()
Reset to starting point of line search.
void moveMesh()
Calculates surface mesh movement.
virtual labelList getActiveDesignVariables() const
Return active design variables.
vectorField controlPointMovement(const scalarField &correction)
List< vectorField > cpsInit_
Backup of initial control points. Useful for line-search.
TypeName("volumetricBSplines")
Runtime type information.
virtual ~optMeshMovementVolumetricBSplines()=default
Destructor.
virtual scalar computeEta(const scalarField &correction)
Compute eta value based on max displacement.
volBSplinesBase & volBSplinesBase_
Reference to underlaying volumetric B-Splines morpher.
Abstract base class for translating an update of the design variables into mesh movement.
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
dynamicFvMesh & mesh
Namespace for OpenFOAM.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73