cyclicTransform.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "cyclicTransform.H"
27 #include "unitConversion.H"
28 #include "IOmanip.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  template<>
36  {
37  "unspecified",
38  "none",
39  "rotational",
40  "translational"
41  };
42 
45 
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::cyclicTransform::update()
53 {
54  if (!transformComplete_)
55  {
56  return;
57  }
58 
59  switch (transformType_)
60  {
61  case UNSPECIFIED:
62  break;
63 
64  case NONE:
65  transform_ = transformer();
66  break;
67 
68  case ROTATIONAL:
69  if (rotationAngle_ == 0)
70  {
71  transform_ = transformer();
72  }
73  else
74  {
75  const tensor R =
76  quaternion(rotationAxis_, degToRad(rotationAngle_)).R();
77 
78  if (mag(rotationCentre_) == 0)
79  {
80  transform_ = transformer::rotation(R);
81  }
82  else
83  {
84  transform_ =
85  transformer::translation(rotationCentre_)
87  & transformer::translation(- rotationCentre_);
88  }
89  }
90  break;
91 
92  case TRANSLATIONAL:
93  if (mag(separation_) == 0)
94  {
95  transform_ = transformer();
96  }
97  else
98  {
99  transform_ = transformer::translation(separation_);
100  }
101  break;
102  }
103 }
104 
105 
106 bool Foam::cyclicTransform::set
107 (
108  const cyclicTransform& t,
109  const scalar lengthScale,
110  const scalar matchTolerance
111 )
112 {
113  // If the supplied transform is unspecified then there is nothing to do
114  if (t.transformType_ == UNSPECIFIED)
115  {
116  return true;
117  }
118 
119  // If this transform is specified then we need to check that it is
120  // compatible with the supplied transform
121  if (transformType_ != UNSPECIFIED)
122  {
123  // If the transforms are of different types then they are incompatible
124  if (transformType_ != t.transformType_)
125  {
126  return false;
127  }
128 
129  // Length-tolerance
130  const scalar lengthTolerance = lengthScale*matchTolerance;
131 
132  // If the transforms are both rotational then the axes must be in the
133  // same direction, the centre points must lie on the same line, and the
134  // angles (if available) must be opposing.
135  if (transformType_ == ROTATIONAL)
136  {
137  const scalar dot = rotationAxis_ & t.rotationAxis_;
138 
139  if (mag(dot) < 1 - matchTolerance)
140  {
141  return false;
142  }
143 
144  if
145  (
146  (rotationAxis_ & (rotationCentre_ - t.rotationCentre_))
147  > lengthTolerance
148  )
149  {
150  return false;
151  }
152 
153  if (transformComplete_ && t.transformComplete_)
154  {
155  if
156  (
157  mag(degToRad(rotationAngle_ - sign(dot)*t.rotationAngle_))
158  > matchTolerance
159  )
160  {
161  return false;
162  }
163  }
164  }
165 
166  // If the transforms are both translational then the separations must
167  // be opposing
168  if (transformType_ == TRANSLATIONAL)
169  {
170  if (transformComplete_ && t.transformComplete_)
171  {
172  if (mag(separation_ - t.separation_) > lengthTolerance)
173  {
174  return false;
175  }
176  }
177  }
178  }
179 
180  // If the supplied transform is more complete then overwrite this with it
181  if (!transformComplete_ && t.transformComplete_)
182  {
183  *this = t;
184  }
185 
186  return true;
187 }
188 
189 
190 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
191 
193 (
194  const bool defaultIsNone
195 )
196 :
197  transformType_(defaultIsNone ? NONE : UNSPECIFIED),
198  rotationAxis_(vector::uniform(NaN)),
199  rotationCentre_(vector::uniform(NaN)),
200  rotationAngle_(NaN),
201  separation_(vector::uniform(NaN)),
202  transformComplete_(transformType_ == NONE),
203  transform_()
204 {}
205 
206 
208 (
209  const dictionary& dict,
210  const bool defaultIsNone
211 )
212 :
213  transformType_
214  (
215  transformTypeNames
216  [
218  (
219  {"transformType", "transform"},
220  transformTypeNames[defaultIsNone ? NONE : UNSPECIFIED]
221  )
222  ]
223  ),
224  rotationAxis_
225  (
226  transformType_ == ROTATIONAL
227  ? normalised(dict.lookup<vector>("rotationAxis"))
228  : vector::uniform(NaN)
229  ),
230  rotationCentre_
231  (
232  transformType_ == ROTATIONAL
233  ? dict.lookup<point>("rotationCentre")
234  : point::uniform(NaN)
235  ),
236  rotationAngle_(dict.lookupOrDefault<scalar>("rotationAngle", NaN)),
237  separation_
238  (
239  transformType_ == TRANSLATIONAL
240  ? (
242  (
243  {"separation", "separationVector"}
244  )
245  )
246  : vector::uniform(NaN)
247  ),
248  transformComplete_
249  (
250  transformType_ == NONE
251  || (
252  transformType_ == ROTATIONAL
253  && dict.found("rotationAngle")
254  )
255  || (
256  transformType_ == TRANSLATIONAL
257  && (dict.found("separation") || dict.found("separationVector"))
258  )
259  ),
260  transform_()
261 {
262  if (transformComplete_)
263  {
264  update();
265  }
266 }
267 
268 
270 (
271  const word& name,
272  const vectorField& areas,
273  const cyclicTransform& transform,
274  const word& nbrName,
275  const cyclicTransform& nbrTransform,
276  const scalar matchTolerance
277 )
278 :
279  cyclicTransform(transform)
280 {
281  // Calculate the total (vector) areas for the supplied patch data
282  const vector area = sum(areas);
283 
284  // Calculate patch length scales
285  const scalar lengthScale = sqrt(mag(area));
286 
287  // Copy as much data from the neighbour as possible
288  if (!transformComplete_ && nbrTransform.transformType_ != UNSPECIFIED)
289  {
290  if (!set(inv(nbrTransform), lengthScale, matchTolerance))
291  {
293  << "Patch " << name
294  << " and neighbour patch " << nbrName
295  << " have incompatible transforms:" << nl << nl << incrIndent;
296 
298  << indent << name << nl << indent << token::BEGIN_BLOCK << nl
299  << incrIndent;
300 
302 
304  << decrIndent << indent << token::END_BLOCK << nl << nl;
305 
307  << indent << nbrName << nl << indent << token::BEGIN_BLOCK << nl
308  << incrIndent;
309 
310  nbrTransform.write(FatalError);
311 
313  << decrIndent << indent << token::END_BLOCK << nl << nl;
314 
316  << decrIndent << exit(FatalError);
317  }
318  }
319 }
320 
321 
323 (
324  const word& name,
325  const pointField& ctrs,
326  const vectorField& areas,
327  const cyclicTransform& transform,
328  const word& nbrName,
329  const pointField& nbrCtrs,
330  const vectorField& nbrAreas,
331  const cyclicTransform& nbrTransform,
332  const scalar matchTolerance
333 )
334 :
336  (
337  name,
338  areas,
339  transform,
340  nbrName,
341  nbrTransform,
342  matchTolerance
343  )
344 {
345  // If there is no geometry from which to calculate the transform then
346  // nothing can be calculated
347  if (areas.size() == 0 || nbrAreas.size() == 0) return;
348 
349  // Calculate the total (vector) areas for the supplied patch data
350  const vector area = sum(areas);
351  const vector nbrArea = sum(nbrAreas);
352 
353  // Calculate the centroids for the supplied patch data
354  const scalarField magAreas(mag(areas));
355  const scalarField magNbrAreas(mag(nbrAreas));
356  const point ctr = sum(ctrs*magAreas)/sum(magAreas);
357  const point nbrCtr = sum(nbrCtrs*magNbrAreas)/sum(magNbrAreas);
358 
359  // Calculate patch length scales
360  const scalar lengthScale = sqrt(sum(magAreas));
361 
362  // Calculate the transformation from the patch geometry
363  if (!transformComplete_)
364  {
365  // Store the old transformation type
366  const transformTypes oldTransformType = transformType_;
367 
368  // Calculate the average patch normals
369  const vector normal = normalised(area);
370  const vector negNbrNormal = - normalised(nbrArea);
371 
372  // Calculate the angle and distance separations
373  const scalar dot = normal & negNbrNormal;
374  const vector delta = ctr - nbrCtr;
375 
376  // Determine the type of transformation if it has not been specified
377  if (transformType_ == UNSPECIFIED)
378  {
379  transformType_ =
380  dot < 1 - rootSmall
381  ? ROTATIONAL
382  : mag(delta) > lengthScale*rootSmall
383  ? TRANSLATIONAL
384  : NONE;
385  }
386 
387  // If the transformation is known to be rotational, then we need to
388  // calculate the angle. If the transformation was previously
389  // unspecified then we also need to calculate the axis and the centre
390  // of rotation.
391  if (transformType_ == ROTATIONAL)
392  {
393  // Calculate the axis, if necessary
394  if (transformType_ != oldTransformType)
395  {
396  const vector midNormal = normalised(normal + negNbrNormal);
397  const vector axis =
398  (ctr - nbrCtr)
399  ^ (
400  normal*(negNbrNormal & midNormal)
401  - negNbrNormal*(normal & midNormal)
402  );
403  const vector axis180 =
404  (ctr - nbrCtr)
405  ^ (normal - negNbrNormal);
406 
407  rotationAxis_ =
408  normalised
409  (
410  mag(axis) > lengthScale*rootSmall
411  ? axis
412  : axis180
413  );
414  }
415 
416  // Calculate the angle
417  const tensor PerpA = tensor::I - sqr(rotationAxis_);
418  const vector normalPerpA = normalised(PerpA & normal);
419  const vector negNbrNormalPerpA = normalised(PerpA & negNbrNormal);
420  const scalar theta =
421  acos(min(max(normalPerpA & negNbrNormalPerpA, -1), 1));
422  rotationAngle_ =
423  - sign((normalPerpA ^ negNbrNormalPerpA) & rotationAxis_)
424  *radToDeg(theta);
425 
426  // Calculate the centre of rotation, if necessary
427  if (transformType_ != oldTransformType)
428  {
429  const tensor R = quaternion(rotationAxis_, theta).R();
430  tensor A = tensor::I - R;
431  vector b = ctr - (R & nbrCtr);
432  const label i = findMax(cmptMag(rotationAxis_));
433  forAll(b, j)
434  {
435  A(i, j) = j == i;
436  }
437  b[i] = 0;
438  rotationCentre_ = inv(A) & b;
439  }
440  }
441 
442  // If the transformation is known to be translational then we just need
443  // to set the separation.
444  if (transformType_ == TRANSLATIONAL)
445  {
446  separation_ = delta;
447  }
448 
449  // Update the transform object
450  transformComplete_ = true;
451  update();
452 
453  // Print results of calculation
454  if (debug)
455  {
456  Info<< "Transformation calculated between patches " << name
457  << " and " << nbrName << ":" << nl << token::BEGIN_BLOCK << nl
458  << incrIndent;
459 
461 
463  }
464  }
465 
466  // Check the transformation is correct to within the matching tolerance
467  const point nbrCtrT =
468  transform_.transformPosition(nbrCtr);
469 
470  const scalar ctrNbrCtrTDistance = mag(ctr - nbrCtrT);
471 
472  if (ctrNbrCtrTDistance > lengthScale*matchTolerance)
473  {
475  << "The distance between the centre of patch " << name
476  << " and the transformed centre of patch " << nbrName << " is "
477  << ctrNbrCtrTDistance << "."
478  << nl
479  << "This is greater than the match tolerance of "
480  << lengthScale*matchTolerance << " for the patch."
481  << nl
482  << "Check that the patches are geometrically similar and that any "
483  << "transformations defined between them are correct"
484  << nl
485  << "It might be possible to fix this problem by increasing the "
486  << "\"matchTolerance\" setting for this patch in the boundary "
487  << exit(FatalError);
488  }
489 }
490 
491 
492 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
493 
495 {}
496 
497 
498 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
499 
501 {
502  const label oldPrecision = os.precision();
503 
504  os.precision(16);
505 
506  if (transformType_ != UNSPECIFIED)
507  {
508  writeEntry(os, "transformType", transformTypeNames[transformType_]);
509  }
510 
511  if (transformType_ == ROTATIONAL)
512  {
513  writeEntry(os, "rotationAxis", rotationAxis_);
514  writeEntry(os, "rotationCentre", rotationCentre_);
515 
516  if (transformComplete_)
517  {
518  writeEntry(os, "rotationAngle", rotationAngle_);
519  }
520  }
521 
522  if (transformType_ == TRANSLATIONAL)
523  {
524  if (transformComplete_)
525  {
526  writeEntry(os, "separation", separation_);
527  }
528  }
529 
530  os.precision(oldPrecision);
531 }
532 
533 
534 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
535 
536 Foam::cyclicTransform Foam::operator&
537 (
538  const transformer& t,
539  const cyclicTransform& c0
540 )
541 {
542  cyclicTransform c1(c0);
543 
544  if (c1.transformType_ == cyclicTransform::ROTATIONAL)
545  {
546  c1.rotationAxis_ = normalised(t.transform(c1.rotationAxis_));
547  c1.rotationCentre_ = t.transformPosition(c1.rotationCentre_);
548  }
549 
550  if (c1.transformType_ == cyclicTransform::TRANSLATIONAL)
551  {
552  if (c1.transformComplete_)
553  {
554  c1.separation_ = t.transform(c1.separation_);
555  }
556  }
557 
558  if (c1.transformComplete_)
559  {
560  c1.update();
561  }
562 
563  return c1;
564 }
565 
566 
568 {
569  cyclicTransform c1(c0);
570 
571  if (c1.transformType_ == cyclicTransform::ROTATIONAL)
572  {
573  if (c1.transformComplete_)
574  {
575  c1.rotationAngle_ = - c1.rotationAngle_;
576  }
577  }
578 
579  if (c1.transformType_ == cyclicTransform::TRANSLATIONAL)
580  {
581  if (c1.transformComplete_)
582  {
583  c1.separation_ = - c1.separation_;
584  }
585  }
586 
587  if (c1.transformComplete_)
588  {
589  c1.update();
590  }
591 
592  return c1;
593 }
594 
595 
596 // ************************************************************************* //
scalar delta
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tensor R() const
The rotation tensor corresponding the quaternion.
Definition: quaternionI.H:323
static transformer translation(const vector &t)
Construct a pure translation transformer.
Definition: transformerI.H:31
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static transformer rotation(const tensor &T)
Construct a pure rotation transformer.
Definition: transformerI.H:43
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
cyclicTransform(const bool defaultIsNone)
void write(Ostream &os) const
Write the data to a dictionary.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
T lookupOrDefaultBackwardsCompatible(const wordList &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, trying a list of keywords in sequence.
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
A class for handling words, derived from string.
Definition: word.H:59
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:60
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
static const NamedEnum< transformTypes, 4 > transformTypeNames
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Istream and Ostream manipulators taking arguments.
static const Tensor I
Definition: Tensor.H:83
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define R(A, B, C, D, E, F, K, M)
messageStream Info
static Vector< scalar > uniform(const scalar &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
dimensioned< scalar > mag(const dimensioned< Type > &)
Cyclic plane tranformation.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:855
Namespace for OpenFOAM.
virtual int precision() const =0
Get precision of output field.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844