globalIndexAndTransformI.H
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) 2011-2018 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 = gi_.processor(a);
37  label procB = gi_.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 = gi_.index(a);
51  label indexB = gi_.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 = gi_.transformIndex(a);
65  label transformB = gi_.transformIndex(b);
66 
67  return transformA < transformB;
68  }
69  }
70 }
71 
72 
74 (
75  const labelList& 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 label transformIndex
112 ) const
113 {
114  labelList permutation(transforms_.size(), 0);
115 
116  label t = transformIndex;
117  forAll(permutation, i)
118  {
119  permutation[i] = (t%3)-1;
120  t /= 3;
121  }
122 
123  return permutation;
124 }
125 
126 
128 (
129  const label transformIndex,
130  const label patchi,
131  const bool isSendingSide,
132  const scalar tol
133 ) const
134 {
135  const labelPair& transSign = patchTransformSign_[patchi];
136 
137  label matchTransI = transSign.first();
138 
139  if (matchTransI >= transforms_.size())
140  {
142  << "patch:" << mesh_.boundaryMesh()[patchi].name()
143  << " transform:" << matchTransI
144  << " out of possible transforms:" << transforms_
145  << exit(FatalError);
146  return labelMin;
147  }
148  else if (matchTransI == -1)
149  {
150  // No additional transformation for this patch
151  return transformIndex;
152  }
153  else
154  {
155  // Decode current set of transforms
156  labelList permutation(decodeTransformIndex(transformIndex));
157 
158 
159  // Add patch transform
160  // ~~~~~~~~~~~~~~~~~~~
161 
162  label sign = transSign.second();
163  if (!isSendingSide)
164  {
165  sign = -sign;
166  }
167 
168 
169  // If this transform been found already by a patch?
170  if (permutation[matchTransI] != 0)
171  {
172  if (sign == 0)
173  {
174  // sent from patch without a transformation. Do nothing.
176  << "patch:" << mesh_.boundaryMesh()[patchi].name()
177  << " transform:" << matchTransI << " sign:" << sign
178  << " current transforms:" << permutation
179  << exit(FatalError);
180  }
181  else if (sign == permutation[matchTransI])
182  {
183  // This is usually illegal. The only exception is for points
184  // on the axis of a 180 degree cyclic wedge when the
185  // transformation is going to be (-1 0 0 0 -1 0 0 0 +1)
186  // (or a different permutation but always two times -1 and
187  // once +1)
188  bool antiCyclic = false;
189 
190  const vectorTensorTransform& vt = transforms_[matchTransI];
191  if (mag(vt.t()) < small && vt.hasR())
192  {
193  const tensor& R = vt.R();
194  scalar sumDiag = tr(R);
195  scalar sumMagDiag = mag(R.xx())+mag(R.yy())+mag(R.zz());
196 
197  if (mag(sumMagDiag-3) < tol && mag(sumDiag+1) < tol)
198  {
199  antiCyclic = true;
200  }
201  }
202 
203  if (antiCyclic)
204  {
205  // 180 degree rotational. Reset transformation.
206  permutation[matchTransI] = 0;
207  }
208  else
209  {
211  << "More than one patch accessing the same transform "
212  << "but not of the same sign." << endl
213  << "patch:" << mesh_.boundaryMesh()[patchi].name()
214  << " transform:" << matchTransI << " sign:" << sign
215  << " current transforms:" << permutation
216  << exit(FatalError);
217  }
218  }
219  else
220  {
221  permutation[matchTransI] = 0;
222  }
223  }
224  else
225  {
226  permutation[matchTransI] = sign;
227  }
228 
229 
230  // Re-encode permutation
231  // ~~~~~~~~~~~~~~~~~~~~~
232 
233  return encodeTransformIndex(permutation);
234  }
235 }
236 
237 
239 (
240  const label transformIndex0,
241  const label transformIndex1
242 ) const
243 {
244  if (transformIndex0 == transformIndex1)
245  {
246  return transformIndex0;
247  }
248 
249 
250  // Count number of transforms
251  labelList permutation0(decodeTransformIndex(transformIndex0));
252  label n0 = 0;
253  forAll(permutation0, i)
254  {
255  if (permutation0[i] != 0)
256  {
257  n0++;
258  }
259  }
260 
261  labelList permutation1(decodeTransformIndex(transformIndex1));
262  label n1 = 0;
263  forAll(permutation1, i)
264  {
265  if (permutation1[i] != 0)
266  {
267  n1++;
268  }
269  }
270 
271  if (n0 <= n1)
272  {
273  return transformIndex0;
274  }
275  else
276  {
277  return transformIndex1;
278  }
279 }
280 
281 
283 (
284  const label transformIndex0,
285  const label transformIndex1
286 ) const
287 {
288  labelList permutation0(decodeTransformIndex(transformIndex0));
289  labelList permutation1(decodeTransformIndex(transformIndex1));
290 
291  forAll(permutation0, i)
292  {
293  permutation0[i] -= permutation1[i];
294  }
295 
296  return encodeTransformIndex(permutation0);
297 }
298 
299 
301 (
302  const label index,
303  const label transformIndex
304 ) const
305 {
306  return encode(Pstream::myProcNo(), index, transformIndex);
307 }
308 
309 
311 (
312  const label proci,
313  const label index,
314  const label transformIndex
315 ) const
316 {
317  if (transformIndex < 0 || transformIndex >= transformPermutations_.size())
318  {
320  << "TransformIndex " << transformIndex
321  << " is outside allowed range of 0 to "
322  << transformPermutations_.size() - 1
323  << abort(FatalError);
324  }
325 
326  if (proci > labelMax/transformPermutations_.size())
327  {
329  << "Overflow : encoding processor " << proci
330  << " in base " << transformPermutations_.size()
331  << " exceeds capability of label (" << labelMax
332  << "). Please recompile with larger datatype for label."
333  << exit(FatalError);
334  }
335 
336  return labelPair
337  (
338  index,
339  transformIndex + proci*transformPermutations_.size()
340  );
341 }
342 
343 
345 (
346  const labelPair& globalIAndTransform
347 ) const
348 {
349  return globalIAndTransform.first();
350 }
351 
352 
354 (
355  const labelPair& globalIAndTransform
356 ) const
357 {
358  return globalIAndTransform.second()/transformPermutations_.size();
359 }
360 
361 
363 (
364  const labelPair& globalIAndTransform
365 ) const
366 {
367  return globalIAndTransform.second()%transformPermutations_.size();
368 }
369 
370 
372 {
373  return transforms_.size();
374 }
375 
376 
379 {
380  return transforms_;
381 }
382 
383 
386 {
387  return transformPermutations_;
388 }
389 
390 
392 {
393  return nullTransformIndex_;
394 }
395 
396 
397 const Foam::labelPairList&
399 {
400  return patchTransformSign_;
401 }
402 
403 
405 (
406  label transformIndex
407 ) const
408 {
409  return transformPermutations_[transformIndex];
410 }
411 
412 
414 (
415  const labelHashSet& patchis
416 ) const
417 {
418  labelList permutation(transforms_.size(), 0);
419 
420  labelList selectedTransformIs(0);
421 
422  if (patchis.empty() || transforms_.empty())
423  {
424  return selectedTransformIs;
425  }
426 
427  forAllConstIter(labelHashSet, patchis, iter)
428  {
429  label patchi = iter.key();
430 
431  const labelPair& transSign = patchTransformSign_[patchi];
432 
433  label matchTransI = transSign.first();
434 
435  if (matchTransI > -1)
436  {
437  label sign = transSign.second();
438 
439  // If this transform been found already by a patch?
440  if (permutation[matchTransI] != 0)
441  {
442  // If so, if they have opposite signs, then this is
443  // considered an error. They are allowed to be the
444  // same sign, but this only results in a single
445  // transform.
446  if (permutation[matchTransI] != sign)
447  {
449  << "More than one patch accessing the same transform "
450  << "but not of the same sign."
451  << exit(FatalError);
452  }
453  }
454  else
455  {
456  permutation[matchTransI] = sign;
457  }
458  }
459  }
460 
461  label nUsedTrans = round(sum(mag(permutation)));
462 
463  if (nUsedTrans == 0)
464  {
465  return selectedTransformIs;
466  }
467 
468  // Number of selected transformations
469  label nSelTrans = pow(label(2), nUsedTrans) - 1;
470 
471  // Pout<< nl << permutation << nl << endl;
472 
473  selectedTransformIs.setSize(nSelTrans);
474 
475  switch (nUsedTrans)
476  {
477  case 1:
478  {
479  selectedTransformIs[0] = encodeTransformIndex(permutation);
480 
481  break;
482  }
483  case 2:
484  {
485  labelList tempPermutation = permutation;
486 
487  label a = 0;
488  label b = 1;
489 
490  // When there are two selected transforms out of three, we
491  // need to choose which of them are being permuted
492  if (transforms_.size() > nUsedTrans)
493  {
494  if (permutation[0] == 0)
495  {
496  a = 1;
497  b = 2;
498  }
499  else if (permutation[1] == 0)
500  {
501  a = 0;
502  b = 2;
503  }
504  else if (permutation[2] == 0)
505  {
506  a = 0;
507  b = 1;
508  }
509  }
510 
511  tempPermutation[a] = a;
512  tempPermutation[b] = permutation[b];
513 
514  selectedTransformIs[0] = encodeTransformIndex(tempPermutation);
515 
516  tempPermutation[a] = permutation[a];
517  tempPermutation[b] = a;
518 
519  selectedTransformIs[1] = encodeTransformIndex(tempPermutation);
520 
521  tempPermutation[a] = permutation[a];
522  tempPermutation[b] = permutation[b];
523 
524  selectedTransformIs[2] = encodeTransformIndex(tempPermutation);
525 
526  break;
527  }
528  case 3:
529  {
530  labelList tempPermutation = permutation;
531 
532  tempPermutation[0] = 0;
533  tempPermutation[1] = 0;
534  tempPermutation[2] = permutation[2];
535 
536  selectedTransformIs[0] = encodeTransformIndex(tempPermutation);
537 
538  tempPermutation[0] = 0;
539  tempPermutation[1] = permutation[1];
540  tempPermutation[2] = 0;
541 
542  selectedTransformIs[1] = encodeTransformIndex(tempPermutation);
543 
544  tempPermutation[0] = 0;
545  tempPermutation[1] = permutation[1];
546  tempPermutation[2] = permutation[2];
547 
548  selectedTransformIs[2] = encodeTransformIndex(tempPermutation);
549 
550  tempPermutation[0] = permutation[0];
551  tempPermutation[1] = 0;
552  tempPermutation[2] = 0;
553 
554  selectedTransformIs[3] = encodeTransformIndex(tempPermutation);
555 
556  tempPermutation[0] = permutation[0];
557  tempPermutation[1] = 0;
558  tempPermutation[2] = permutation[2];
559 
560  selectedTransformIs[4] = encodeTransformIndex(tempPermutation);
561 
562  tempPermutation[0] = permutation[0];
563  tempPermutation[1] = permutation[1];
564  tempPermutation[2] = 0;
565 
566  selectedTransformIs[5] = encodeTransformIndex(tempPermutation);
567 
568  tempPermutation[0] = permutation[0];
569  tempPermutation[1] = permutation[1];
570  tempPermutation[2] = permutation[2];
571 
572  selectedTransformIs[6] = encodeTransformIndex(tempPermutation);
573 
574  break;
575  }
576  default:
577  {
579  << "Only 1-3 transforms are possible."
580  << exit(FatalError);
581  }
582  }
583 
584  return selectedTransformIs;
585 }
586 
587 
589 (
590  const labelHashSet& patchis,
591  const point& pt
592 ) const
593 {
594  labelList transIs = transformIndicesForPatches(patchis);
595 
596  // Pout<< patchis << nl << transIs << endl;
597 
598  pointField transPts(transIs.size());
599 
600  forAll(transIs, tII)
601  {
602  transPts[tII] = transformPermutations_[transIs[tII]].transformPosition
603  (
604  pt
605  );
606  }
607 
608  return transPts;
609 }
610 
611 
612 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
label index(const labelPair &globalIAndTransform) const
Index carried by the object.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label nIndependentTransforms() const
Return the number of independent transforms.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Vector-tensor class used to perform translations and rotations in 3D space.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
const Cmpt & yy() const
Definition: TensorI.H:181
const Cmpt & xx() const
Definition: TensorI.H:153
const List< vectorTensorTransform > & transformPermutations() const
Return access to the permuted transforms.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
label processor(const labelPair &globalIAndTransform) const
Which processor does this come from?
labelPair encode(const label index, const label transformIndex) const
Encode index and bare index as components on own processor.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
label addToTransformIndex(const label transformIndex, const label patchi, const bool isSendingSide=true, const scalar tol=small) const
Add patch transformation to transformIndex. Return new.
label subtractTransformIndex(const label transformIndex0, const label transformIndex1) const
Subtract two transformIndices.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const List< vectorTensorTransform > & transforms() const
Return access to the stored independent transforms.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
label nullTransformIndex() const
Return the transformIndex (index in transformPermutations)
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
pointField transformPatches(const labelHashSet &patchIs, const point &pt) const
Apply all of the transform permutations.
static const label labelMax
Definition: label.H:62
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const Type & second() const
Return second.
Definition: Pair.H:99
label transformIndex(const labelPair &globalIAndTransform) const
Transform carried by the object.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
#define R(A, B, C, D, E, F, K, M)
const vectorTensorTransform & transform(label transformIndex) const
Access the overall (permuted) transform corresponding.
const List< labelPair > & patchTransformSign() const
Return access to the per-patch transform-sign pairs.
label minimumTransformIndex(const label transformIndex0, const label transformIndex1) const
Combine two transformIndices.
dimensioned< scalar > mag(const dimensioned< Type > &)
labelList transformIndicesForPatches(const labelHashSet &patchIs) const
Access the all of the indices of the transform.
const Type & first() const
Return first.
Definition: Pair.H:87
const Cmpt & zz() const
Definition: TensorI.H:209
labelList decodeTransformIndex(const label transformIndex) const
Decode transform index.
static const label labelMin
Definition: label.H:61
label encodeTransformIndex(const labelList &permutationIndices) const
Generate a transform index from the permutation indices of.