noConstraint.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) 2021 PCOpt/NTUA
9  Copyright (C) 2021 FOSS GP
10 -------------------------------------------------------------------------------
11 License
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 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "noConstraint.H"
32 
33 // * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(noConstraint, 1);
39  (
40  morphingBoxConstraint,
41  noConstraint,
42  dictionary
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 (
51  scalarField& dvSens,
52  const scalarField& cpSens
53 )
54 {
55  dvSens = cpSens;
56 }
57 
58 
60 (
61  autoPtr<scalarField>& lowerBounds,
62  autoPtr<scalarField>& upperBounds,
63  const NURBS3DVolume& boxI,
64  const label passed
65 )
66 {
67  const vectorField& cps = boxI.getControlPoints();
68  const Vector<label> nCPsDir = boxI.nCPsPerDirection();
69  // Internal points
70  for (label k = 1; k < nCPsDir[2] - 1; ++k)
71  {
72  for (label j = 1; j < nCPsDir[1] - 1; ++j)
73  {
74  for (label i = 1; i < nCPsDir[0] - 1; ++i)
75  {
76  label cpID(boxI.getCPID(i, j, k));
77  for (label idir = 0; idir < 3; ++idir)
78  {
79  label iIncr(idir == 0);
80  label jIncr(idir == 1);
81  label kIncr(idir == 2);
82  label prevCP
83  (boxI.getCPID(i - iIncr, j - jIncr, k - kIncr));
84  label nextCP
85  (boxI.getCPID(i + iIncr, j + jIncr, k + kIncr));
86  lowerBounds()[3*cpID + idir + passed] =
87  0.5
88  *(
89  cps[prevCP].component(idir)
90  + cps[cpID].component(idir)
91  );
92  upperBounds()[3*cpID + idir + passed] =
93  0.5
94  *(
95  cps[nextCP].component(idir)
96  + cps[cpID].component(idir)
97  );
98  }
99  }
100  }
101  }
102 }
103 
104 
106 (
107  autoPtr<scalarField>& lowerBounds,
108  autoPtr<scalarField>& upperBounds,
109  const NURBS3DVolume& boxI,
110  const label passed
111 )
112 {
113  const vectorField& cps = boxI.getControlPoints();
114  const Vector<label> nCPsDir = boxI.nCPsPerDirection();
115  // Loop over boundaries in all directions of the box
116  for (label ibound = 0; ibound < 3; ++ibound)
117  {
118  // Start of iterators in the three directions
119  Vector<label> minID(1, 1, 1);
120  // End of iterators in the three directions
121  Vector<label> maxID(nCPsDir[0] - 2, nCPsDir[1] - 2, nCPsDir[2] - 2);
122  // Increment of iterators in the three directions
123  Vector<label> incr(1, 1, 1);
124 
125  // Adjust looping in the direction we are examining
126  minID[ibound] = 0;
127  maxID[ibound] = nCPsDir[ibound];
128  incr[ibound] = nCPsDir[ibound] - 1;
129  Vector<label> indices(Zero);
130  label& i = indices[0];
131  label& j = indices[1];
132  label& k = indices[2];
133 
134  for (k = minID[2]; k < maxID[2]; k += incr[2])
135  {
136  for (j = minID[1]; j < maxID[1]; j += incr[1])
137  {
138  for (i = minID[0]; i < maxID[0]; i += incr[0])
139  {
140  label cpID(boxI.getCPID(i, j, k));
141  for (label dir = 0; dir < 3; ++dir)
142  {
143  Vector<label> incrMinus(dir == 0, dir == 1, dir == 2);
144  Vector<label> incrPlus(dir == 0, dir == 1, dir == 2);
145  // Adjust increment for the ibound direction
146  incrMinus[ibound] =
147  label
148  (
149  incrMinus[ibound]
150  && indices[ibound] == nCPsDir[ibound] - 1
151  );
152  incrPlus[ibound] =
153  label(incrMinus[ibound] && indices[ibound] == 0);
154  label prevCP =
155  boxI.getCPID
156  (
157  i - incrMinus[0],
158  j - incrMinus[1],
159  k - incrMinus[2]
160  );
161  label nextCP =
162  boxI.getCPID
163  (
164  i + incrPlus[0],
165  j + incrPlus[1],
166  k + incrPlus[2]
167  );
168  if (incrMinus[ibound])
169  {
170  lowerBounds()[3*cpID + dir + passed] =
171  0.5
172  *(
173  cps[prevCP].component(dir)
174  + cps[cpID].component(dir)
175  );
176  }
177  if (incrPlus[ibound])
178  {
179  upperBounds()[3*cpID + dir + passed] =
180  0.5
181  *(
182  cps[nextCP].component(dir)
183  + cps[cpID].component(dir)
184  );
185  }
186  }
187  }
188  }
189  }
190  }
191 }
192 
193 
194 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
195 
196 Foam::noConstraint::noConstraint
197 (
198  const fvMesh& mesh,
199  const dictionary& dict,
200  volumetricBSplinesDesignVariables& designVariables
201 )
202 :
203  morphingBoxConstraint(mesh, dict, designVariables)
204 {
206 }
207 
208 
209 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
210 
212 (
213  autoPtr<scalarField>& lowerBounds,
214  autoPtr<scalarField>& upperBounds
215 )
216 {
217  // Does nothing
218 }
219 
220 
222 (
223  autoPtr<scalarField>& lowerBounds,
224  autoPtr<scalarField>& upperBounds
225 )
226 {
227  if (designVariables_.nonOverlappingCPs() && designVariables_.updateBounds())
228  {
229  DebugInfo
230  << "Updating bounds for the design variables " << endl;
231  const PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxesRef();
232  label passed(0);
233  for (const NURBS3DVolume& boxI : boxes)
234  {
235  // Bounds for internal control points
236  updateInternalBounds(lowerBounds, upperBounds, boxI, passed);
237  // Bounds for boundary points.
238  // Assumes that the boundary edges remain fixed
239  updateBoundaryBounds(lowerBounds, upperBounds, boxI, passed);
240  passed += 3*boxI.getControlPoints().size();
241  }
242  DebugInfo
243  << "lower bounds " << lowerBounds() << endl;
245  << "upper bounds " << upperBounds() << endl;
246  }
247 }
248 
249 
251 (
252  const labelList& activeCPCoors
253 )
254 {
255  return activeCPCoors;
256 }
257 
258 
260 (
262 )
263 {
264  return designVariables;
265 }
266 
267 
269 (
270  const scalarField& cps
271 )
272 {
273  return cps;
274 }
275 
276 
278 (
279  const scalarField& correction
280 )
281 {
282  return correction;
283 }
284 
285 
286 // ************************************************************************* //
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
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual tmp< scalarField > correctionCPs(const scalarField &correctionDVs)
Convert the correction of the design variables to the correction of the control points.
Definition: noConstraint.C:271
virtual labelList computeActiveDesignVariables(const labelList &activeCPCoors)
Compute the active design variables based on the IDs of the active control point coordinates.
Definition: noConstraint.C:244
virtual void updateBounds(autoPtr< scalarField > &lowerBounds, autoPtr< scalarField > &upperBounds)
Update the bounds of the design variables.
Definition: noConstraint.C:215
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
void updateBoundaryBounds(autoPtr< scalarField > &lowerBounds, autoPtr< scalarField > &upperBounds, const NURBS3DVolume &boxI, const label passed)
Update the bounds of the boundary control points.
Definition: noConstraint.C:99
label k
Boltzmann constant.
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:69
virtual tmp< scalarField > controlPointsToDesignVariables(const scalarField &cps)
Return the design variables corresponding to the given control points.
Definition: noConstraint.C:262
volBSplinesBase & volBSplinesBase_
Easy access to the volBSplinesBase resting in the designVariables_.
Macros for easy insertion into run-time selection tables.
virtual void computeDVsSensitivities(scalarField &dvSens, const scalarField &cpSens)
Compute sensitivities wrt the design variables (chain rule)
Definition: noConstraint.C:43
Vector< label > nCPsPerDirection() const
Get number of control points per direction.
void updateInternalBounds(autoPtr< scalarField > &lowerBounds, autoPtr< scalarField > &upperBounds, const NURBS3DVolume &boxI, const label passed)
Update the bounds of the internal control points.
Definition: noConstraint.C:53
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
virtual tmp< scalarField > designVariablesToControlPoints(const scalarField &designVariables)
Convert design variables to control points, stored in a scalarField.
Definition: noConstraint.C:253
label getCPID(const label i, const label j, const label k) const
Get control point ID from its I-J-K coordinates.
#define DebugInfo
Report an information message using Foam::Info.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:608
defineTypeNameAndDebug(combustionModel, 0)
virtual void computeBounds(autoPtr< scalarField > &lowerBounds, autoPtr< scalarField > &upperBounds)
Transform bounds from control points to design variables.
Definition: noConstraint.C:205
volumetricBSplinesDesignVariables & designVariables_
Reference to underlaying volumetric B-Splines morpher.
label getTotalControlPointsNumber() const
Get cumulative number of control points from all boxes.
Abstract base class for defining design variables.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Field< vector > vectorField
Specialisation of Field<T> for vector.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const vectorField & getControlPoints() const
Get control points.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127