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-2022 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  cyclicTransform(true)
195 {}
196 
197 
199 (
200  const bool defaultIsNone
201 )
202 :
203  transformType_(defaultIsNone ? NONE : UNSPECIFIED),
204  rotationAxis_(vector::uniform(NaN)),
205  rotationCentre_(vector::uniform(NaN)),
206  rotationAngle_(NaN),
207  separation_(vector::uniform(NaN)),
208  transformComplete_(transformType_ == NONE),
209  transform_()
210 {}
211 
212 
214 (
215  const dictionary& dict,
216  const bool defaultIsNone
217 )
218 :
219  transformType_
220  (
222  [
224  (
225  {"transformType", "transform"},
226  transformTypeNames[defaultIsNone ? NONE : UNSPECIFIED]
227  )
228  ]
229  ),
230  rotationAxis_
231  (
232  transformType_ == ROTATIONAL
233  ? normalised(dict.lookup<vector>("rotationAxis"))
234  : vector::uniform(NaN)
235  ),
236  rotationCentre_
237  (
238  transformType_ == ROTATIONAL
239  ? dict.lookup<point>("rotationCentre")
240  : point::uniform(NaN)
241  ),
242  rotationAngle_(dict.lookupOrDefault<scalar>("rotationAngle", NaN)),
243  separation_
244  (
245  transformType_ == TRANSLATIONAL
246  ? (
248  (
249  {"separation", "separationVector"}
250  )
251  )
252  : vector::uniform(NaN)
253  ),
254  transformComplete_
255  (
256  transformType_ == NONE
257  || (
258  transformType_ == ROTATIONAL
259  && dict.found("rotationAngle")
260  )
261  || (
262  transformType_ == TRANSLATIONAL
263  && (dict.found("separation") || dict.found("separationVector"))
264  )
265  ),
266  transform_()
267 {
268  if (transformComplete_)
269  {
270  update();
271  }
272 }
273 
274 
276 (
277  const word& name,
278  const vectorField& areas,
279  const cyclicTransform& transform,
280  const word& nbrName,
281  const cyclicTransform& nbrTransform,
282  const scalar matchTolerance
283 )
284 :
285  cyclicTransform(transform)
286 {
287  // Calculate the total (vector) areas for the supplied patch data
288  const vector area = sum(areas);
289 
290  // Calculate patch length scales
291  const scalar lengthScale = sqrt(mag(area));
292 
293  // Copy as much data from the neighbour as possible
294  if (!transformComplete_ && nbrTransform.transformType_ != UNSPECIFIED)
295  {
296  if (!set(inv(nbrTransform), lengthScale, matchTolerance))
297  {
299  << "Patch " << name
300  << " and neighbour patch " << nbrName
301  << " have incompatible transforms:" << nl << nl << incrIndent;
302 
304  << indent << name << nl << indent << token::BEGIN_BLOCK << nl
305  << incrIndent;
306 
308 
310  << decrIndent << indent << token::END_BLOCK << nl << nl;
311 
313  << indent << nbrName << nl << indent << token::BEGIN_BLOCK << nl
314  << incrIndent;
315 
316  nbrTransform.write(FatalError);
317 
319  << decrIndent << indent << token::END_BLOCK << nl << nl;
320 
322  << decrIndent << exit(FatalError);
323  }
324  }
325 }
326 
327 
329 (
330  const word& name,
331  const pointField& ctrs,
332  const vectorField& areas,
333  const cyclicTransform& transform,
334  const word& nbrName,
335  const pointField& nbrCtrs,
336  const vectorField& nbrAreas,
337  const cyclicTransform& nbrTransform,
338  const scalar matchTolerance
339 )
340 :
342  (
343  name,
344  areas,
345  transform,
346  nbrName,
347  nbrTransform,
348  matchTolerance
349  )
350 {
351  // If there is no geometry from which to calculate the transform then
352  // nothing can be calculated
353  if (areas.size() == 0 || nbrAreas.size() == 0) return;
354 
355  // Calculate the total (vector) areas for the supplied patch data
356  const vector area = sum(areas);
357  const vector nbrArea = sum(nbrAreas);
358 
359  // Calculate the centroids for the supplied patch data
360  const scalarField magAreas(mag(areas));
361  const scalarField magNbrAreas(mag(nbrAreas));
362  const point ctr = sum(ctrs*magAreas)/sum(magAreas);
363  const point nbrCtr = sum(nbrCtrs*magNbrAreas)/sum(magNbrAreas);
364 
365  // Calculate patch length scales
366  const scalar lengthScale = sqrt(sum(magAreas));
367 
368  // Calculate the transformation from the patch geometry
369  if (!transformComplete_)
370  {
371  // Store the old transformation type
372  const transformTypes oldTransformType = transformType_;
373 
374  // Calculate the average patch normals
375  const vector normal = normalised(area);
376  const vector negNbrNormal = - normalised(nbrArea);
377 
378  // Calculate the angle and distance separations
379  const scalar dot = normal & negNbrNormal;
380  const vector delta = ctr - nbrCtr;
381 
382  // Determine the type of transformation if it has not been specified
383  if (transformType_ == UNSPECIFIED)
384  {
385  transformType_ =
386  dot < 1 - rootSmall
387  ? ROTATIONAL
388  : mag(delta) > lengthScale*rootSmall
389  ? TRANSLATIONAL
390  : NONE;
391  }
392 
393  // If the transformation is known to be rotational, then we need to
394  // calculate the angle. If the transformation was previously
395  // unspecified then we also need to calculate the axis and the centre
396  // of rotation.
397  if (transformType_ == ROTATIONAL)
398  {
399  // Calculate the axis, if necessary
400  if (transformType_ != oldTransformType)
401  {
402  const vector midNormal = normalised(normal + negNbrNormal);
403  const vector axis =
404  (ctr - nbrCtr)
405  ^ (
406  normal*(negNbrNormal & midNormal)
407  - negNbrNormal*(normal & midNormal)
408  );
409  const vector axis180 =
410  (ctr - nbrCtr)
411  ^ (normal - negNbrNormal);
412 
413  rotationAxis_ =
414  normalised
415  (
416  mag(axis) > lengthScale*rootSmall
417  ? axis
418  : axis180
419  );
420  }
421 
422  // Calculate the angle
423  const tensor PerpA = tensor::I - sqr(rotationAxis_);
424  const vector normalPerpA = normalised(PerpA & normal);
425  const vector negNbrNormalPerpA = normalised(PerpA & negNbrNormal);
426  const scalar theta =
427  acos(min(max(normalPerpA & negNbrNormalPerpA, -1), 1));
428  rotationAngle_ =
429  - sign((normalPerpA ^ negNbrNormalPerpA) & rotationAxis_)
430  *radToDeg(theta);
431 
432  // Calculate the centre of rotation, if necessary
433  if (transformType_ != oldTransformType)
434  {
435  const tensor R = quaternion(rotationAxis_, theta).R();
436  tensor A = tensor::I - R;
437  vector b = ctr - (R & nbrCtr);
438  const label i = findMax(cmptMag(rotationAxis_));
439  forAll(b, j)
440  {
441  A(i, j) = j == i;
442  }
443  b[i] = 0;
444  rotationCentre_ = inv(A) & b;
445  }
446  }
447 
448  // If the transformation is known to be translational then we just need
449  // to set the separation.
450  if (transformType_ == TRANSLATIONAL)
451  {
452  separation_ = delta;
453  }
454 
455  // Update the transform object
456  transformComplete_ = true;
457  update();
458 
459  // Print results of calculation
460  if (debug)
461  {
462  Info<< "Transformation calculated between patches " << name
463  << " and " << nbrName << ":" << nl << token::BEGIN_BLOCK << nl
464  << incrIndent;
465 
467 
469  }
470  }
471 
472  // Check the transformation is correct to within the matching tolerance
473  const point nbrCtrT =
474  transform_.transformPosition(nbrCtr);
475 
476  const scalar ctrNbrCtrTDistance = mag(ctr - nbrCtrT);
477 
478  if (ctrNbrCtrTDistance > lengthScale*matchTolerance)
479  {
481  << "The distance between the centre of patch " << name
482  << " and the transformed centre of patch " << nbrName << " is "
483  << ctrNbrCtrTDistance << "."
484  << nl
485  << "This is greater than the match tolerance of "
486  << lengthScale*matchTolerance << " for the patch."
487  << nl
488  << "Check that the patches are geometrically similar and that any "
489  << "transformations defined between them are correct"
490  << nl
491  << "It might be possible to fix this problem by increasing the "
492  << "\"matchTolerance\" setting for this patch in the boundary "
493  << exit(FatalError);
494  }
495 }
496 
497 
498 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
499 
501 {}
502 
503 
504 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
505 
507 {
508  const label oldPrecision = os.precision();
509 
510  os.precision(16);
511 
512  if (transformType_ != UNSPECIFIED)
513  {
514  writeEntry(os, "transformType", transformTypeNames[transformType_]);
515  }
516 
517  if (transformType_ == ROTATIONAL)
518  {
519  writeEntry(os, "rotationAxis", rotationAxis_);
520  writeEntry(os, "rotationCentre", rotationCentre_);
521 
522  if (transformComplete_)
523  {
524  writeEntry(os, "rotationAngle", rotationAngle_);
525  }
526  }
527 
528  if (transformType_ == TRANSLATIONAL)
529  {
530  if (transformComplete_)
531  {
532  writeEntry(os, "separation", separation_);
533  }
534  }
535 
536  os.precision(oldPrecision);
537 }
538 
539 
540 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
541 
542 Foam::cyclicTransform Foam::operator&
543 (
544  const transformer& t,
545  const cyclicTransform& c0
546 )
547 {
548  cyclicTransform c1(c0);
549 
550  if (c1.transformType_ == cyclicTransform::ROTATIONAL)
551  {
552  c1.rotationAxis_ = normalised(t.transform(c1.rotationAxis_));
553  c1.rotationCentre_ = t.transformPosition(c1.rotationCentre_);
554  }
555 
556  if (c1.transformType_ == cyclicTransform::TRANSLATIONAL)
557  {
558  if (c1.transformComplete_)
559  {
560  c1.separation_ = t.transform(c1.separation_);
561  }
562  }
563 
564  if (c1.transformComplete_)
565  {
566  c1.update();
567  }
568 
569  return c1;
570 }
571 
572 
574 {
575  cyclicTransform c1(c0);
576 
577  if (c1.transformType_ == cyclicTransform::ROTATIONAL)
578  {
579  if (c1.transformComplete_)
580  {
581  c1.rotationAngle_ = - c1.rotationAngle_;
582  }
583  }
584 
585  if (c1.transformType_ == cyclicTransform::TRANSLATIONAL)
586  {
587  if (c1.transformComplete_)
588  {
589  c1.separation_ = - c1.separation_;
590  }
591  }
592 
593  if (c1.transformComplete_)
594  {
595  c1.update();
596  }
597 
598  return c1;
599 }
600 
601 
602 // ************************************************************************* //
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:663
#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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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
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
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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).
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
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 > &)
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
Cyclic plane transformation.
friend cyclicTransform inv(const cyclicTransform &c)
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:875
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:864