globalIndexAndTransformI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 "polyMesh.H"
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 bool Foam::globalIndexAndTransform::less::operator()
31 (
32  const labelPair& a,
33  const labelPair& b
34 ) const
35 {
36  label procA = globalIndexAndTransform::processor(a);
37  label procB = globalIndexAndTransform::processor(b);
38 
39  if (procA < procB)
40  {
41  return true;
42  }
43  else if (procA > procB)
44  {
45  return false;
46  }
47  else
48  {
49  // Equal proc.
50  label indexA = globalIndexAndTransform::index(a);
51  label indexB = globalIndexAndTransform::index(b);
52 
53  if (indexA < indexB)
54  {
55  return true;
56  }
57  else if (indexA > indexB)
58  {
59  return false;
60  }
61  else
62  {
63  // Equal index
64  label transformA = globalIndexAndTransform::transformIndex(a);
65  label transformB = globalIndexAndTransform::transformIndex(b);
66 
67  return transformA < transformB;
68  }
69  }
70 }
71 
72 
74 (
75  const List<label>& permutationIndices
76 ) const
77 {
78  if (permutationIndices.size() != transforms_.size())
79  {
81  << "permutationIndices " << permutationIndices
82  << "are of a different size to the number of independent transforms"
83  << abort(FatalError);
84  }
85 
86  label transformIndex = 0;
87 
88  label w = 1;
89 
90  forAll(transforms_, b)
91  {
92  if (mag(permutationIndices[b]) > 1)
93  {
95  << "permutationIndices " << permutationIndices
96  << "are illegal, they must all be only -1, 0 or +1"
97  << abort(FatalError);
98  }
99 
100  transformIndex += (permutationIndices[b] + 1)*w;
101 
102  w *= 3;
103  }
104 
105  return transformIndex;
106 }
107 
108 
110 (
111  const FixedList<Foam::label, 3>& permutation
112 ) const
113 {
114  if (nIndependentTransforms() == 0)
115  {
116  return 0;
117  }
118  if (nIndependentTransforms() == 1)
119  {
120  return permutation[0]+1;
121  }
122  else if (nIndependentTransforms() == 2)
123  {
124  return (permutation[1]+1)*3 + (permutation[0]+1);
125  }
126  else
127  {
128  return
129  (permutation[2]+1)*9
130  + (permutation[1]+1)*3
131  + (permutation[0]+1);
132  }
133 }
134 
135 
138 (
139  const label transformIndex
140 ) const
141 {
142  FixedList<label, 3> permutation(label(0));
143 
144  label t = transformIndex;
145  if (nIndependentTransforms() > 0)
146  {
147  permutation[0] = (t%3)-1;
148  if (nIndependentTransforms() > 1)
149  {
150  t /= 3;
151  permutation[1] = (t%3)-1;
152  if (nIndependentTransforms() > 2)
153  {
154  t /= 3;
155  permutation[2] = (t%3)-1;
156  }
157  }
158  }
159 
160 # ifdef FULLDEBUG
161  t /= 3;
162  if (t != 0)
163  {
165  << "transformIndex : " << transformIndex
166  << " has more than 3 fields."
167  << abort(FatalError);
168  }
169 # endif
170 
171  return permutation;
172 }
173 
174 
176 (
177  const label transformIndex,
178  const label patchI,
179  const bool isSendingSide,
180  const scalar tol
181 ) const
182 {
183  const Pair<label>& transSign = patchTransformSign_[patchI];
184 
185  label matchTransI = transSign.first();
186 
187  // Hardcoded for max 3 transforms only!
188 
189  if (matchTransI > -1 && matchTransI < 3)
190  {
191  FixedList<label, 3> permutation = decodeTransformIndex(transformIndex);
192 
193 
194  // Add patch transform
195  // ~~~~~~~~~~~~~~~~~~~
196 
197  label sign = transSign.second();
198  if (!isSendingSide)
199  {
200  sign = -sign;
201  }
202 
203 
204  // If this transform been found already by a patch?
205  if (permutation[matchTransI] != 0)
206  {
207  if (sign == 0)
208  {
209  // sent from patch without a transformation. Do nothing.
211  << "patch:" << mesh_.boundaryMesh()[patchI].name()
212  << " transform:" << matchTransI << " sign:" << sign
213  << " current transforms:" << permutation
214  << exit(FatalError);
215  }
216  else if (sign == permutation[matchTransI])
217  {
218  // This is usually illegal. The only exception is for points
219  // on the axis of a 180 degree cyclic wedge when the
220  // transformation is going to be (-1 0 0 0 -1 0 0 0 +1)
221  // (or a different permutation but always two times -1 and
222  // once +1)
223  bool antiCyclic = false;
224 
225  const vectorTensorTransform& vt = transforms_[matchTransI];
226  if (mag(vt.t()) < SMALL && vt.hasR())
227  {
228  const tensor& R = vt.R();
229  scalar sumDiag = tr(R);
230  scalar sumMagDiag = mag(R.xx())+mag(R.yy())+mag(R.zz());
231 
232  if (mag(sumMagDiag-3) < tol && mag(sumDiag+1) < tol)
233  {
234  antiCyclic = true;
235  }
236  }
237 
238  if (antiCyclic)
239  {
240  // 180 degree rotational. Reset transformation.
241  permutation[matchTransI] = 0;
242  }
243  else
244  {
246  << "More than one patch accessing the same transform "
247  << "but not of the same sign." << endl
248  << "patch:" << mesh_.boundaryMesh()[patchI].name()
249  << " transform:" << matchTransI << " sign:" << sign
250  << " current transforms:" << permutation
251  << exit(FatalError);
252  }
253  }
254  else
255  {
256  permutation[matchTransI] = 0;
257  }
258  }
259  else
260  {
261  permutation[matchTransI] = sign;
262  }
263 
264 
265  // Re-encode permutation
266  // ~~~~~~~~~~~~~~~~~~~~~
267 
268  return encodeTransformIndex(permutation);
269  }
270  else
271  {
272  return transformIndex;
273  }
274 }
275 
276 
278 (
279  const label transformIndex0,
280  const label transformIndex1
281 ) const
282 {
283  if (transformIndex0 == transformIndex1)
284  {
285  return transformIndex0;
286  }
287 
288 
289  // Count number of transforms
290  FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
291  label n0 = 0;
292  forAll(permutation0, i)
293  {
294  if (permutation0[i] != 0)
295  {
296  n0++;
297  }
298  }
299 
300  FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);
301  label n1 = 0;
302  forAll(permutation1, i)
303  {
304  if (permutation1[i] != 0)
305  {
306  n1++;
307  }
308  }
309 
310  if (n0 <= n1)
311  {
312  return transformIndex0;
313  }
314  else
315  {
316  return transformIndex1;
317  }
318 }
319 
320 
322 (
323  const label transformIndex0,
324  const label transformIndex1
325 ) const
326 {
327  FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
328  FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);
329 
330  forAll(permutation0, i)
331  {
332  permutation0[i] -= permutation1[i];
333  }
334 
335  return encodeTransformIndex(permutation0);
336 }
337 
338 
340 (
341  const label index,
342  const label transformIndex
343 )
344 {
345  return encode(Pstream::myProcNo(), index, transformIndex);
346 }
347 
348 
350 (
351  const label procI,
352  const label index,
353  const label transformIndex
354 )
355 {
356  if (transformIndex < 0 || transformIndex >= base_)
357  {
359  << "TransformIndex " << transformIndex
360  << " is outside allowed range of 0 to "
361  << base_ - 1
362  << abort(FatalError);
363  }
364 
365  if (procI > labelMax/base_)
366  {
368  << "Overflow : encoding processor " << procI << " in base " << base_
369  << " exceeds capability of label (" << labelMax
370  << "). Please recompile with larger datatype for label."
371  << exit(FatalError);
372  }
373 
374  return labelPair
375  (
376  index,
377  transformIndex + procI*base_
378  );
379 }
380 
381 
383 (
384  const labelPair& globalIAndTransform
385 )
386 {
387  return globalIAndTransform.first();
388 }
389 
390 
392 (
393  const labelPair& globalIAndTransform
394 )
395 {
396  return globalIAndTransform.second()/base_;
397 }
398 
399 
401 (
402  const labelPair& globalIAndTransform
403 )
404 {
405  return globalIAndTransform.second() % base_;
406 }
407 
408 
410 {
411  return transforms_.size();
412 }
413 
414 
417 {
418  return transforms_;
419 }
420 
421 
424 {
425  return transformPermutations_;
426 }
427 
428 
430 {
431  return nullTransformIndex_;
432 }
433 
434 
437 {
438  return patchTransformSign_;
439 }
440 
441 
443 (
444  label transformIndex
445 ) const
446 {
447  return transformPermutations_[transformIndex];
448 }
449 
450 
452 (
453  const labelHashSet& patchIs
454 ) const
455 {
456  List<label> permutation(transforms_.size(), 0);
457 
458  labelList selectedTransformIs(0);
459 
460  if (patchIs.empty() || transforms_.empty())
461  {
462  return selectedTransformIs;
463  }
464 
465  forAllConstIter(labelHashSet, patchIs, iter)
466  {
467  label patchI = iter.key();
468 
469  const Pair<label>& transSign = patchTransformSign_[patchI];
470 
471  label matchTransI = transSign.first();
472 
473  if (matchTransI > -1)
474  {
475  label sign = transSign.second();
476 
477  // If this transform been found already by a patch?
478  if (permutation[matchTransI] != 0)
479  {
480  // If so, if they have opposite signs, then this is
481  // considered an error. They are allowed to be the
482  // same sign, but this only results in a single
483  // transform.
484  if (permutation[matchTransI] != sign)
485  {
487  << "More than one patch accessing the same transform "
488  << "but not of the same sign."
489  << exit(FatalError);
490  }
491  }
492  else
493  {
494  permutation[matchTransI] = sign;
495  }
496  }
497  }
498 
499  label nUsedTrans = round(sum(mag(permutation)));
500 
501  if (nUsedTrans == 0)
502  {
503  return selectedTransformIs;
504  }
505 
506  // Number of selected transformations
507  label nSelTrans = pow(label(2), nUsedTrans) - 1;
508 
509  // Pout<< nl << permutation << nl << endl;
510 
511  selectedTransformIs.setSize(nSelTrans);
512 
513  switch (nUsedTrans)
514  {
515  case 1:
516  {
517  selectedTransformIs[0] = encodeTransformIndex(permutation);
518 
519  break;
520  }
521  case 2:
522  {
523  List<label> tempPermutation = permutation;
524 
525  label a = 0;
526  label b = 1;
527 
528  // When there are two selected transforms out of three, we
529  // need to choose which of them are being permuted
530  if (transforms_.size() > nUsedTrans)
531  {
532  if (permutation[0] == 0)
533  {
534  a = 1;
535  b = 2;
536  }
537  else if (permutation[1] == 0)
538  {
539  a = 0;
540  b = 2;
541  }
542  else if (permutation[2] == 0)
543  {
544  a = 0;
545  b = 1;
546  }
547  }
548 
549  tempPermutation[a] = a;
550  tempPermutation[b] = permutation[b];
551 
552  selectedTransformIs[0] = encodeTransformIndex(tempPermutation);
553 
554  tempPermutation[a] = permutation[a];
555  tempPermutation[b] = a;
556 
557  selectedTransformIs[1] = encodeTransformIndex(tempPermutation);
558 
559  tempPermutation[a] = permutation[a];
560  tempPermutation[b] = permutation[b];
561 
562  selectedTransformIs[2] = encodeTransformIndex(tempPermutation);
563 
564  break;
565  }
566  case 3:
567  {
568  List<label> tempPermutation = permutation;
569 
570  tempPermutation[0] = 0;
571  tempPermutation[1] = 0;
572  tempPermutation[2] = permutation[2];
573 
574  selectedTransformIs[0] = encodeTransformIndex(tempPermutation);
575 
576  tempPermutation[0] = 0;
577  tempPermutation[1] = permutation[1];
578  tempPermutation[2] = 0;
579 
580  selectedTransformIs[1] = encodeTransformIndex(tempPermutation);
581 
582  tempPermutation[0] = 0;
583  tempPermutation[1] = permutation[1];
584  tempPermutation[2] = permutation[2];
585 
586  selectedTransformIs[2] = encodeTransformIndex(tempPermutation);
587 
588  tempPermutation[0] = permutation[0];
589  tempPermutation[1] = 0;
590  tempPermutation[2] = 0;
591 
592  selectedTransformIs[3] = encodeTransformIndex(tempPermutation);
593 
594  tempPermutation[0] = permutation[0];
595  tempPermutation[1] = 0;
596  tempPermutation[2] = permutation[2];
597 
598  selectedTransformIs[4] = encodeTransformIndex(tempPermutation);
599 
600  tempPermutation[0] = permutation[0];
601  tempPermutation[1] = permutation[1];
602  tempPermutation[2] = 0;
603 
604  selectedTransformIs[5] = encodeTransformIndex(tempPermutation);
605 
606  tempPermutation[0] = permutation[0];
607  tempPermutation[1] = permutation[1];
608  tempPermutation[2] = permutation[2];
609 
610  selectedTransformIs[6] = encodeTransformIndex(tempPermutation);
611 
612  break;
613  }
614  default:
615  {
617  << "Only 1-3 transforms are possible."
618  << exit(FatalError);
619  }
620  }
621 
622  return selectedTransformIs;
623 }
624 
625 
627 (
628  const labelHashSet& patchIs,
629  const point& pt
630 ) const
631 {
632  labelList transIs = transformIndicesForPatches(patchIs);
633 
634  // Pout<< patchIs << nl << transIs << endl;
635 
636  pointField transPts(transIs.size());
637 
638  forAll(transIs, tII)
639  {
640  transPts[tII] = transformPermutations_[transIs[tII]].transformPosition
641  (
642  pt
643  );
644  }
645 
646  return transPts;
647 }
648 
649 
650 // ************************************************************************* //
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::globalIndexAndTransform::decodeTransformIndex
FixedList< label, 3 > decodeTransformIndex(const label transformIndex) const
Decode transform index. Hardcoded to 3 independent transforms max.
Definition: globalIndexAndTransformI.H:138
Foam::vectorTensorTransform::R
const tensor & R() const
Definition: vectorTensorTransformI.H:73
Foam::globalIndexAndTransform::transformPermutations
const List< vectorTensorTransform > & transformPermutations() const
Return access to the permuted transforms.
Definition: globalIndexAndTransformI.H:423
Foam::labelMax
static const label labelMax
Definition: label.H:62
Foam::globalIndexAndTransform::nIndependentTransforms
label nIndependentTransforms() const
Return the number of independent transforms.
Definition: globalIndexAndTransformI.H:409
Foam::Pair::first
const Type & first() const
Return first.
Definition: Pair.H:87
Foam::globalIndexAndTransform::encode
static labelPair encode(const label index, const label transformIndex)
Encode index and bare index as components on own processor.
Definition: globalIndexAndTransformI.H:340
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:179
polyMesh.H
Foam::HashSet< label, Hash< label > >
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::globalIndexAndTransform::patchTransformSign
const List< Pair< label > > & patchTransformSign() const
Return access to the per-patch transform-sign pairs.
Definition: globalIndexAndTransformI.H:436
R
#define R(A, B, C, D, E, F, K, M)
Foam::globalIndexAndTransform::nullTransformIndex
label nullTransformIndex() const
Return the transformIndex (index in transformPermutations)
Definition: globalIndexAndTransformI.H:429
Foam::label
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
Foam::globalIndexAndTransform::transformIndex
static label transformIndex(const labelPair &globalIAndTransform)
Transform carried by the object.
Definition: globalIndexAndTransformI.H:401
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::HashTable::empty
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
Foam::globalIndexAndTransform::index
static label index(const labelPair &globalIAndTransform)
Index carried by the object.
Definition: globalIndexAndTransformI.H:383
Foam::Pair::second
const Type & second() const
Return second.
Definition: Pair.H:99
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::vectorTensorTransform::t
const vector & t() const
Definition: vectorTensorTransformI.H:67
Foam::FatalError
error FatalError
Foam::globalIndexAndTransform::transformIndicesForPatches
labelList transformIndicesForPatches(const labelHashSet &patchIs) const
Access the all of the indices of the transform.
Definition: globalIndexAndTransformI.H:452
Foam::globalIndexAndTransform::encodeTransformIndex
label encodeTransformIndex(const FixedList< Foam::label, 3 > &permutationIndices) const
Encode transform index. Hardcoded to 3 independent transforms max.
Definition: globalIndexAndTransformI.H:110
Foam::globalIndexAndTransform::transformPatches
pointField transformPatches(const labelHashSet &patchIs, const point &pt) const
Apply all of the transform permutations.
Definition: globalIndexAndTransformI.H:627
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::vectorTensorTransform::hasR
bool hasR() const
Definition: vectorTensorTransformI.H:79
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::globalIndexAndTransform::subtractTransformIndex
label subtractTransformIndex(const label transformIndex0, const label transformIndex1) const
Subtract two transformIndices.
Definition: globalIndexAndTransformI.H:322
Foam::globalIndexAndTransform::transform
const vectorTensorTransform & transform(label transformIndex) const
Access the overall (permuted) transform corresponding.
Definition: globalIndexAndTransformI.H:443
Foam::List::setSize
void setSize(const label)
Reset size of List.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::globalIndexAndTransform::processor
static label processor(const labelPair &globalIAndTransform)
Which processor does this come from?
Definition: globalIndexAndTransformI.H:392
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
Foam::Vector< scalar >
Foam::List< label >
Foam::vectorTensorTransform
Vector-tensor class used to perform translations and rotations in 3D space.
Definition: vectorTensorTransform.H:61
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::globalIndexAndTransform::transforms_
List< vectorTensorTransform > transforms_
The possible independent (non-permuted) transforms of the.
Definition: globalIndexAndTransform.H:93
Foam::globalIndexAndTransform::addToTransformIndex
label addToTransformIndex(const label transformIndex, const label patchI, const bool isSendingSide=true, const scalar tol=SMALL) const
Add patch transformation to transformIndex. Return new.
Definition: globalIndexAndTransformI.H:176
Foam::globalIndexAndTransform::transforms
const List< vectorTensorTransform > & transforms() const
Return access to the stored independent transforms.
Definition: globalIndexAndTransformI.H:416
Foam::IOstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
Foam::globalIndexAndTransform::minimumTransformIndex
label minimumTransformIndex(const label transformIndex0, const label transformIndex1) const
Combine two transformIndices.
Definition: globalIndexAndTransformI.H:278
Foam::labelPair
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48