Loading...
Searching...
No Matches
volumetricBSplinesDesignVariables.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-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::volumetricBSplinesDesignVariables
29
30Description
31 Volumetric B-Splines design variables for shape optimisation
32
33SourceFiles
34 volumetricBSplinesDesignVariables.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef volumetricBSplinesDesignVariables_H
39#define volumetricBSplinesDesignVariables_H
40
42#include "volBSplinesBase.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
50/*---------------------------------------------------------------------------*\
51 Class volumetricBSplinesDesignVariables Declaration
52\*---------------------------------------------------------------------------*/
53
54class volumetricBSplinesDesignVariables
55:
56 public shapeDesignVariables,
58{
59protected:
60
61 // Protected data
62
63 //- Reference to underlaying volumetric B-Splines morpher
65
66 //- Should the control points be non-overlapping
68
69 //- Should the bounds of the non-overlapping control points be updated
70 //- in each optimisation cycle
71 bool updateBounds_;
73 //- Constraint imposed on the movement of the design variables.
74 // Can be used to e.g. impose a uniform movement of the control points
75 // in one direction, etc.
77
78
79 // Protected Member Functions
81 //- Size of the active control points
82 virtual label sensSize() const;
83
84 //- Components of the active control points
85 virtual const labelList& activeSensitivities() const;
86
87 //- Set IDs of active design variables
88 // Might be different than what volBSplinesBase_ returns, if
89 // constraint_ changes the number of design variables
91
92 //- Set control points based on current design variables values
94
95 //- Set the design variables based on the current control points
97
98 //- Set the design variables based on the given control points
99 void controlPointsToDesignVariables(const vectorField& controlPoints);
100
101 //- Read bounds for design variables, if present
102 void readBounds
103 (
104 autoPtr<scalar> lowerBoundPtr = nullptr,
105 autoPtr<scalar> upperBoundPtr = nullptr
106 );
107
108 //- Read one set of bounds (lower, upper)
109 void readBounds
110 (
112 const word& boundsName,
113 const label sign
114 );
115
116 //- Write current bounds to file
117 void writeBounds(const scalarField& bounds, const word& name) const;
118
119 //- Set the field driving the displacement method.
120 // Can be either the movement of the control points or the boundary
121 // displacement, depending on the method used to move the mesh
122 void setDisplacement(const vectorField& cpMovement);
123
124
125private:
126
127 // Private Member Functions
128
129 //- Disallow default bitwise copy construct
130 volumetricBSplinesDesignVariables
131 (
132 const volumetricBSplinesDesignVariables&
133 ) = delete;
134
135 //- Disallow default bitwise assignment
136 void operator=(const volumetricBSplinesDesignVariables&) = delete;
137
138
139public:
140
141 //- Runtime type information
142 TypeName("volumetricBSplines");
143
144
145 // Constructors
146
147 //- Construct from components
148 volumetricBSplinesDesignVariables
149 (
150 fvMesh& mesh,
151 const dictionary& dict
152 );
153
154
155 //- Destructor
156 virtual ~volumetricBSplinesDesignVariables() = default;
157
158
159 // Member Functions
160
161 //- Are control points non-overlapping
162 inline bool nonOverlappingCPs() const
163 {
164 return nonOverlappingCPs_;
165 }
166
167 //- Are bounds to be updated after each optmisation cycle
168 inline bool updateBounds() const
169 {
170 return updateBounds_;
171 }
172
173 //- Constant access to the volBSplinesBase object
175 {
176 return volBSplinesBase_;
177 }
178
179 //- Non-constant access to the volBSplinesBase object
181 {
182 return volBSplinesBase_;
183 }
184
185 //- Transform correction of design variables to control points movement
186 tmp<vectorField> controlPointMovement
187 (
189 );
190
191 //- Update design variables based on a given correction
193
194 //- Reset to starting point of line search
195 virtual void resetDesignVariables();
196
197 //- Compute eta if not set in the first step
198 virtual scalar computeEta(scalarField& correction);
199
200 //- Whether to use global sum when computing matrix-vector products
201 //- in update methods
202 virtual bool globalSum() const;
203
204 //- Assemble the sensitivity derivatives, by also applying possible
205 //- constraints
207 (
209 );
210
211 //- For design variables with a dynamic character (i.e. changing
212 // number), perform the evolution.
213 // Hijacked here to evolve the bounds of the design variables
214 // after the end of each optimisation cycle
215 virtual void evolveNumber();
217 //- Write fields to support continuation
218 virtual bool writeData(Ostream& os) const;
219
220 // Fields related to sensitivity computations
221
222 //- Get dxdb for all mesh points
225 const label varID
226 ) const;
227
228 //- Get dxdb for given design variable and patch
230 (
231 const label patchI,
232 const label varID
233 ) const;
234
235 //- Get dndb for given design variable and patch
236 virtual tmp<vectorField> dndb
237 (
238 const label patchI,
239 const label varID
240 ) const;
241
242 //- Get dSdb for given design variable and patch
243 virtual tmp<vectorField> dSdb
244 (
245 const label patchI,
246 const label varID
247 ) const;
248
249 //- Get dCdb for given design variable.
250 // Used for FI-based sensitivities
251 virtual tmp<volVectorField> dCdb(const label varID) const;
252};
253
254
255// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256
257} // End namespace Foam
258
259// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260
261#endif
262
263// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Abstract base class for adjoint-based sensitivities.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
localIOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
A class for managing temporary objects.
Definition tmp.H:75
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
virtual void evolveNumber()
For design variables with a dynamic character (i.e. changing.
virtual tmp< vectorField > dSdb(const label patchI, const label varID) const
Get dSdb for given design variable and patch.
virtual ~volumetricBSplinesDesignVariables()=default
Destructor.
virtual scalar computeEta(scalarField &correction)
Compute eta if not set in the first step.
virtual void resetDesignVariables()
Reset to starting point of line search.
const volBSplinesBase & getVolBSplinesBase() const
Constant access to the volBSplinesBase object.
bool updateBounds() const
Are bounds to be updated after each optmisation cycle.
void controlPointsToDesignVariables()
Set the design variables based on the current control points.
void setDisplacement(const vectorField &cpMovement)
Set the field driving the displacement method.
virtual bool writeData(Ostream &os) const
Write fields to support continuation.
void readBounds(autoPtr< scalar > lowerBoundPtr=nullptr, autoPtr< scalar > upperBoundPtr=nullptr)
Read bounds for design variables, if present.
void designVariablesToControlPoints()
Set control points based on current design variables values.
void setActiveDesignVariables()
Set IDs of active design variables.
TypeName("volumetricBSplines")
Runtime type information.
virtual bool globalSum() const
Whether to use global sum when computing matrix-vector products in update methods.
bool nonOverlappingCPs() const
Are control points non-overlapping.
void writeBounds(const scalarField &bounds, const word &name) const
Write current bounds to file.
bool updateBounds_
Should the bounds of the non-overlapping control points be updated in each optimisation cycle.
tmp< vectorField > controlPointMovement(const scalarField &correction)
Transform correction of design variables to control points movement.
volBSplinesBase & getVolBSplinesBase()
Non-constant access to the volBSplinesBase object.
virtual tmp< vectorField > dxdbVol(const label varID) const
Get dxdb for all mesh points.
virtual const labelList & activeSensitivities() const
Components of the active control points.
bool nonOverlappingCPs_
Should the control points be non-overlapping.
volBSplinesBase & volBSplinesBase_
Reference to underlaying volumetric B-Splines morpher.
virtual label sensSize() const
Size of the active control points.
autoPtr< morphingBoxConstraint > constraint_
Constraint imposed on the movement of the design variables.
virtual tmp< vectorField > dxdbFace(const label patchI, const label varID) const
Get dxdb for given design variable and patch.
virtual tmp< vectorField > dndb(const label patchI, const label varID) const
Get dndb for given design variable and patch.
virtual tmp< volVectorField > dCdb(const label varID) const
Get dCdb for given design variable.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Assemble the sensitivity derivatives, by also applying possible constraints.
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh update()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for bounding specifications. At the moment, mostly for tables.
Namespace for OpenFOAM.
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68