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