Bezier.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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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 
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "pointFields.H"
33 #include "Bezier.H"
34 #include "Pstream.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(Bezier, 0);
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 Bezier::Bezier(const fvMesh& mesh, const dictionary& dict)
48 :
49  mesh_(mesh),
50  dict_(dict),
51  nBezier_(dict_.subDict("Bezier").get<label>("nBezier")),
52  dxidXj_(nBezier_),
53  confineXmovement_
54  (
55  dict_.subDict("Bezier").get<boolList>("confineXmovement")
56  ),
57  confineYmovement_
58  (
59  dict_.subDict("Bezier").get<boolList>("confineYmovement")
60  ),
61  confineZmovement_
62  (
63  dict_.subDict("Bezier").get<boolList>("confineZmovement")
64  ),
65  confineMovement_(3, boolList(nBezier_, false)),
66  activeDesignVariables_(3*nBezier_)
67 {
68  if
69  (
70  confineXmovement_.size() != nBezier_
71  || confineYmovement_.size() != nBezier_
72  || confineZmovement_.size() != nBezier_
73  )
74  {
76  << "confineMovement lists sizes "
77  << confineXmovement_.size() << " "
78  << confineYmovement_.size() << " "
79  << confineZmovement_.size() << " "
80  << "are incompatible with nBezier " << nBezier_
81  << endl << endl
82  << exit(FatalError);
83  }
87 
88  // Determine active design variables
89  label iActive(0);
90  for (label iDir = 0; iDir < 3; ++iDir)
91  {
92  for (label iCP = 0; iCP < nBezier_; ++iCP)
93  {
94  if (confineMovement_[iDir][iCP])
95  {
96  activeDesignVariables_[iActive++] = iDir*nBezier_ + iCP;
97  }
98  }
99  }
101 
102  // Read bezier parameterization
103  for (label iCP = 0; iCP < nBezier_; ++iCP)
104  {
105  dxidXj_.set
106  (
107  iCP,
108  new pointTensorField
109  (
110  IOobject
111  (
112  "dxidXj_"+name(iCP),
113  mesh_.time().timeName(),
114  mesh_,
117  ),
119  )
120  );
121  }
122 }
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
127 label Bezier::nBezier() const
128 {
129  return nBezier_;
130 }
131 
132 
134 {
135  return dxidXj_;
136 }
137 
138 
140 {
141  return confineXmovement_;
142 }
143 
144 
146 {
147  return confineYmovement_;
148 }
149 
150 
152 {
153  return confineZmovement_;
154 }
155 
156 
158 {
159  return confineMovement_;
160 }
161 
162 
164 (
165  const label patchI,
166  const label cpI,
167  bool returnDimensionedNormalSens
168 ) const
169 {
170  const fvPatch& patch = mesh_.boundary()[patchI];
171  const polyPatch& ppatch = patch.patch();
172 
173  // Return field
174  auto tdndbSens = tmp<tensorField>::New(patch.size(), Zero);
175  auto& dndbSens = tdndbSens.ref();
176 
177  // Auxiliary quantities
179  const label patchStart = ppatch.start();
180  const tensorField& dxdbInt = dxidXj_[cpI].primitiveField();
181 
182  // Loop over patch faces
183  forAll(patch, fI)
184  {
185  const face& fGlobal = mesh_.faces()[fI + patchStart];
186  const pointField facePoints = fGlobal.points(mesh_.points());
187 
188  // Loop over face points
189  tensorField facePointDerivs(facePoints.size(), Zero);
190  forAll(fGlobal, pI)
191  {
192  facePointDerivs[pI] = dxdbInt[fGlobal[pI]];
193  }
194 
195  // Determine whether to return variance of dimensioned or unit normal
196  if (returnDimensionedNormalSens)
197  {
198  dndbSens[fI] =
200  (
201  facePoints,
202  facePointDerivs
203  )[1];
204  }
205  else
206  {
207  dndbSens[fI] =
209  (
210  facePoints,
211  facePointDerivs
212  )[2];
213  }
214  }
215  return tdndbSens;
216 }
217 
218 
220 (
221  const label patchI,
222  const label cpI,
223  const label idir,
224  bool returnDimensionedNormalSens
225 ) const
226 {
227  const fvPatch& patch = mesh_.boundary()[patchI];
228  const polyPatch& ppatch = patch.patch();
229 
230  // Return field
231  auto tdndbSens = tmp<vectorField>::New(patch.size(), Zero);
232  auto& dndbSens = tdndbSens.ref();
233 
234  // Auxiliary quantities
236  const label patchStart = ppatch.start();
237  const tensorField& dxdbInt = dxidXj_[cpI].primitiveField();
238  vectorField dxdbDir(dxdbInt.size(), Zero);
239  dxdbDir.replace(0, dxdbInt.component(3*idir));
240  dxdbDir.replace(1, dxdbInt.component(3*idir + 1));
241  dxdbDir.replace(2, dxdbInt.component(3*idir + 2));
242 
243  // Loop over patch faces
244  forAll(patch, fI)
245  {
246  const face& fGlobal = mesh_.faces()[fI + patchStart];
247  const pointField facePoints = fGlobal.points(mesh_.points());
248 
249  // Loop over face points
250  vectorField facePointDerivs(facePoints.size(), Zero);
251  forAll(fGlobal, pI)
252  {
253  facePointDerivs[pI] = dxdbDir[fGlobal[pI]];
254  }
255 
256  // Determine whether to return variance of dimensioned or unit normal
257  if (returnDimensionedNormalSens)
258  {
259  dndbSens[fI] =
261  (
262  facePoints, facePointDerivs
263  )[1];
264  }
265  else
266  {
267  dndbSens[fI] =
269  (
270  facePoints,
271  facePointDerivs
272  )[2];
273  }
274  }
275 
276  return tdndbSens;
277 }
278 
279 
281 (
282  const label patchI,
283  const label cpI,
284  bool useChainRule
285 ) const
286 {
287  const polyPatch& patch = mesh_.boundary()[patchI].patch();
288 
289  // Return field
290  auto tdxdbFace = tmp<tensorField>::New(patch.size(), Zero);
291  auto& dxdbFace = tdxdbFace.ref();
292  if (useChainRule)
293  {
294  // Auxiliary quantities
296  const label patchStart = patch.start();
297  const tensorField& dxdbInt = dxidXj_[cpI].primitiveField();
298 
299  // Loop over patch faces
300  forAll(patch, fI)
301  {
302  const face& fGlobal = mesh_.faces()[fI + patchStart];
303  const pointField facePoints = fGlobal.points(mesh_.points());
304 
305  // Loop over face points
306  tensorField facePointDerivs(facePoints.size(), Zero);
307  forAll(fGlobal, pI)
308  {
309  facePointDerivs[pI] = dxdbInt[fGlobal[pI]];
310  }
311  dxdbFace[fI] =
313  (
314  facePoints,
315  facePointDerivs
316  )[0];
317  }
318  }
319  // Simple average of dxdb in points. Older and less accurate
320  else
321  {
323  dxdbFace = patchInter.pointToFaceInterpolate
324  (
325  dxidXj_[cpI].boundaryField()[patchI].patchInternalField()()
326  )();
327  }
328  return tdxdbFace;
329 }
330 
331 
333 (
334  const label patchI,
335  const label cpI,
336  const label idir,
337  bool useChainRule
338 ) const
339 {
340  const polyPatch& patch = mesh_.boundary()[patchI].patch();
341 
342  // Return field
343  auto tdxdbFace = tmp<vectorField>::New(patch.size(), Zero);
344  auto& dxdbFace = tdxdbFace.ref();
345 
346  if (useChainRule)
347  {
348  // Auxiliary quantities
350  const label patchStart = patch.start();
351  const tensorField& dxdbInt = dxidXj_[cpI].primitiveField();
352  vectorField dxdbDir(dxdbInt.size(), Zero);
353  dxdbDir.replace(0, dxdbInt.component(3*idir));
354  dxdbDir.replace(1, dxdbInt.component(3*idir + 1));
355  dxdbDir.replace(2, dxdbInt.component(3*idir + 2));
356 
357  // Loop over patch faces
358  forAll(patch, fI)
359  {
360  const face& fGlobal = mesh_.faces()[fI + patchStart];
361  const pointField facePoints = fGlobal.points(mesh_.points());
362 
363  // Loop over face points
364  vectorField facePointDerivs(facePoints.size(), Zero);
365  forAll(fGlobal, pI)
366  {
367  facePointDerivs[pI] = dxdbDir[fGlobal[pI]];
368  }
369  dxdbFace[fI] =
371  (
372  facePoints,
373  facePointDerivs
374  )[0];
375  }
376  }
377  // Simple average of dxdb in points. Older and less accurate
378  else
379  {
381  vectorField dxdb(patch.nPoints(), Zero);
382  dxdb.replace
383  (
384  idir,
385  dxidXj_[cpI].boundaryField()[patchI].patchInternalField()().
386  component(0)
387  );
388  dxdbFace = patchInter.pointToFaceInterpolate(dxdb)();
389  }
390  return tdxdbFace;
391 }
392 
393 
395 (
396  const label globalFaceI,
397  const label cpI
398 ) const
399 {
400  const face& faceI(mesh_.faces()[globalFaceI]);
401  tensorField fPoints_d(faceI.size(), Zero);
402  forAll(faceI, fpI)
403  {
404  fPoints_d[fpI] = dxidXj_[cpI].primitiveField()[faceI[fpI]];
405  }
406  return fPoints_d;
407 }
408 
409 
411 (
412  const label globalFaceI,
413  const label cpI,
414  const label idir
415 ) const
416 {
417  const face& faceI(mesh_.faces()[globalFaceI]);
418  vectorField fPoints_d(faceI.size(), Zero);
419  forAll(faceI, fpI)
420  {
421  const tensor& dxdbTensor = dxidXj_[cpI].primitiveField()[faceI[fpI]];
422  fPoints_d[fpI].x() = dxdbTensor.component(3*idir);
423  fPoints_d[fpI].y() = dxdbTensor.component(3*idir + 1);
424  fPoints_d[fpI].z() = dxdbTensor.component(3*idir + 2);
425  }
426  return fPoints_d;
427 }
428 
429 
431 {
432  return activeDesignVariables_;
433 }
434 
435 
436 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
437 
438 } // End namespace Foam
439 
440 // ************************************************************************* //
volFields.H
Foam::Tensor< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Bezier::dxidXj_
PtrList< pointTensorField > dxidXj_
Definition: Bezier.H:80
Foam::Bezier::confineMovement_
boolListList confineMovement_
Definition: Bezier.H:85
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::Bezier::confineXmovement_
boolList confineXmovement_
Definition: Bezier.H:82
Foam::deltaBoundary
Differentiation of the mesh data structure.
Definition: deltaBoundary.H:58
Foam::Bezier::confineXmovement
const boolList & confineXmovement() const
Confine x movement.
Definition: Bezier.C:139
Foam::Tensor::z
Vector< Cmpt > z() const
Extract vector for row 2.
Definition: TensorI.H:292
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::PrimitivePatchInterpolation::pointToFaceInterpolate
tmp< Field< Type > > pointToFaceInterpolate(const Field< Type > &pf) const
Interpolate from points to faces.
Definition: PrimitivePatchInterpolation.C:234
Foam::Bezier::getActiveDesignVariables
const labelList & getActiveDesignVariables() const
Return active design variables.
Definition: Bezier.C:430
Foam::Bezier::confineYmovement_
boolList confineYmovement_
Definition: Bezier.H:83
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
Foam::Bezier::confineZmovement
const boolList & confineZmovement() const
Confine z movement.
Definition: Bezier.C:151
Foam::Bezier::activeDesignVariables_
labelList activeDesignVariables_
Definition: Bezier.H:86
Bezier.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Bezier::confineZmovement_
boolList confineZmovement_
Definition: Bezier.H:84
Foam::Bezier::dxidXj
PtrList< pointTensorField > & dxidXj()
dx/db tensor for all control points
Definition: Bezier.C:133
Foam::Bezier::nBezier
label nBezier() const
Number of Bezier control points.
Definition: Bezier.C:127
Foam::Field< tensor >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::Tensor::x
Vector< Cmpt > x() const
Extract vector for row 0.
Definition: TensorI.H:278
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::Field::replace
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:551
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::Bezier::dxdbFace
tmp< tensorField > dxdbFace(const label patchI, const label cpI, bool useChainRule=true) const
dxdb tensor for a Bezier parameterized patch
Definition: Bezier.C:281
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::deltaBoundary::makeFaceCentresAndAreas_d
vectorField makeFaceCentresAndAreas_d(const pointField &p, const pointField &p_d)
Definition: deltaBoundary.C:72
Pstream.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:539
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:361
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::Bezier::facePoints_d
tensorField facePoints_d(const label globalFaceI, const label cpI) const
For a given (global) face ID, return the change of the face points.
Definition: Bezier.C:395
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Tensor::y
Vector< Cmpt > y() const
Extract vector for row 1.
Definition: TensorI.H:285
Foam::PrimitivePatchInterpolation
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
Definition: PrimitivePatchInterpolation.H:53
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< bool >
Foam::Bezier::mesh_
const fvMesh & mesh_
Definition: Bezier.H:76
Foam::face::points
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition: faceI.H:87
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::Bezier::confineYmovement
const boolList & confineYmovement() const
Confine y movement.
Definition: Bezier.C:145
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::Bezier::confineMovement
const boolListList & confineMovement() const
Info about confining movement in all directions.
Definition: Bezier.C:157
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointFields.H
Foam::Bezier::dndbBasedSensitivities
tmp< tensorField > dndbBasedSensitivities(const label patchI, const label cpI, bool returnDimensionedNormalSens=true) const
Definition: Bezier.C:164
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::Bezier::nBezier_
label nBezier_
Definition: Bezier.H:79