BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / OpenFOAM / primitives / globalIndexAndTransform / globalIndexAndTransformI.H
blobc1cdbf91c13ce0ee8a0120e54ee9ea19c98d9271
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "polyMesh.H"
28 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
30 bool Foam::globalIndexAndTransform::less::operator()
32     const labelPair& a,
33     const labelPair& b
34 ) const
36     label procA = globalIndexAndTransform::processor(a);
37     label procB = globalIndexAndTransform::processor(b);
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);
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);
67             return transformA < transformB;
68         }
69     }
73 Foam::label Foam::globalIndexAndTransform::encodeTransformIndex
75     const List<label>& permutationIndices
76 ) const
78     if (permutationIndices.size() != transforms_.size())
79     {
80         FatalErrorIn
81         (
82             "Foam::label encodeTransformIndex"
83             "("
84                 "const List<label>& permutationIndices,"
85             ") const"
86         )
87             << "permutationIndices " << permutationIndices
88             << "are of a different size to the number of independent transforms"
89             << abort(FatalError);
90     }
92     label transformIndex = 0;
94     label w = 1;
96     forAll(transforms_, b)
97     {
98         if (mag(permutationIndices[b]) > 1)
99         {
100             FatalErrorIn
101             (
102                 "Foam::label encodeTransformIndex"
103                 "("
104                 "const List<label>& permutationIndices,"
105                 ") const"
106             )
107                 << "permutationIndices " << permutationIndices
108                 << "are illegal, they must all be only -1, 0 or +1"
109                 << abort(FatalError);
110         }
112         transformIndex += (permutationIndices[b] + 1)*w;
114         w *= 3;
115     }
117     return transformIndex;
121 Foam::label Foam::globalIndexAndTransform::encodeTransformIndex
123     const FixedList<Foam::label, 3>& permutation
124 ) const
126     if (nIndependentTransforms() == 0)
127     {
128         return 0;
129     }
130     if (nIndependentTransforms() == 1)
131     {
132         return permutation[0]+1;
133     }
134     else if (nIndependentTransforms() == 2)
135     {
136         return (permutation[1]+1)*3 + (permutation[0]+1);
137     }
138     else
139     {
140         return
141             (permutation[2]+1)*9
142           + (permutation[1]+1)*3
143           + (permutation[0]+1);
144     }
148 Foam::FixedList<Foam::label, 3>
149 Foam::globalIndexAndTransform::decodeTransformIndex
151     const label transformIndex
152 ) const
154     FixedList<label, 3> permutation(0);
156     label t = transformIndex;
157     if (nIndependentTransforms() > 0)
158     {
159         permutation[0] = (t%3)-1;
160         if (nIndependentTransforms() > 1)
161         {
162             t /= 3;
163             permutation[1] = (t%3)-1;
164             if (nIndependentTransforms() > 2)
165             {
166                 t /= 3;
167                 permutation[2] = (t%3)-1;
168             }
169         }
170     }
172 #   ifdef FULLDEBUG
173     t /= 3;
174     if (t != 0)
175     {
176         FatalErrorIn
177         (
178             "globalIndexAndTransform::decodeTransformIndex(const label)"
179         )   << "transformIndex : " << transformIndex
180             << " has more than 3 fields."
181             << abort(FatalError);
182     }
183 #   endif
185     return permutation;
189 Foam::label Foam::globalIndexAndTransform::addToTransformIndex
191     const label transformIndex,
192     const label patchI,
193     const bool isSendingSide
194 ) const
196     const Pair<label>& transSign = patchTransformSign_[patchI];
198     label matchTransI = transSign.first();
200     // Hardcoded for max 3 transforms only!
202     if (matchTransI > -1 && matchTransI < 3)
203     {
204         FixedList<label, 3> permutation = decodeTransformIndex(transformIndex);
207         // Add patch transform
208         // ~~~~~~~~~~~~~~~~~~~
210         label sign = transSign.second();
211         if (!isSendingSide)
212         {
213             sign = -sign;
214         }
217         // If this transform been found already by a patch?
218         if (permutation[matchTransI] != 0)
219         {
220             if (sign == 0)
221             {
222                 // sent from patch without a transformation. Do nothing.
223                 FatalErrorIn("globalIndexAndTransform::addToTransformIndex(..)")
224                     << "patch:" << mesh_.boundaryMesh()[patchI].name()
225                     << " transform:" << matchTransI << " sign:" << sign
226                     << "  current transforms:" << permutation
227                     << exit(FatalError);
228             }
229             else if (sign == permutation[matchTransI])
230             {
231                 FatalErrorIn
232                 (
233                     "Foam::label "
234                     "Foam::globalIndexAndTransform::addToTransformIndex\n"
235                     "(\n"
236                         "const label,\n"
237                         "const label,\n"
238                         "const bool\n"
239                     ") const\n"
240                 )   << "More than one patch accessing the same transform "
241                     << "but not of the same sign." << endl
242                     << "patch:" << mesh_.boundaryMesh()[patchI].name()
243                     << " transform:" << matchTransI << " sign:" << sign
244                     << "  current transforms:" << permutation
245                     << exit(FatalError);
246             }
247             else
248             {
249                 permutation[matchTransI] = 0;
250             }
251         }
252         else
253         {
254             permutation[matchTransI] = sign;
255         }
258         // Re-encode permutation
259         // ~~~~~~~~~~~~~~~~~~~~~
261         return encodeTransformIndex(permutation);
262     }
263     else
264     {
265         return transformIndex;
266     }
270 Foam::label Foam::globalIndexAndTransform::minimumTransformIndex
272     const label transformIndex0,
273     const label transformIndex1
274 ) const
276     if (transformIndex0 == transformIndex1)
277     {
278         return transformIndex0;
279     }
282     // Count number of transforms
283     FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
284     label n0 = 0;
285     forAll(permutation0, i)
286     {
287         if (permutation0[i] != 0)
288         {
289             n0++;
290         }
291     }
293     FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);
294     label n1 = 0;
295     forAll(permutation1, i)
296     {
297         if (permutation1[i] != 0)
298         {
299             n1++;
300         }
301     }
303     if (n0 <= n1)
304     {
305         return transformIndex0;
306     }
307     else
308     {
309         return transformIndex1;
310     }
314 Foam::label Foam::globalIndexAndTransform::subtractTransformIndex
316     const label transformIndex0,
317     const label transformIndex1
318 ) const
320     FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
321     FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);
323     forAll(permutation0, i)
324     {
325         permutation0[i] -= permutation1[i];
326     }
328     return encodeTransformIndex(permutation0);
332 Foam::labelPair Foam::globalIndexAndTransform::encode
334     const label index,
335     const label transformIndex
338     return encode(Pstream::myProcNo(), index, transformIndex);
342 Foam::labelPair Foam::globalIndexAndTransform::encode
344     const label procI,
345     const label index,
346     const label transformIndex
349     if (transformIndex < 0 || transformIndex >= base_)
350     {
351         FatalErrorIn
352         (
353             "Foam::labelPair Foam::globalIndexAndTransform::encode"
354             "("
355                 "const label procI, "
356                 "const label index, "
357                 "const label transformIndex"
358             ")"
359         )
360             << "TransformIndex " << transformIndex
361             << " is outside allowed range of 0 to "
362             << base_ - 1
363             << abort(FatalError);
364     }
366     if (procI > labelMax/base_)
367     {
368         FatalErrorIn
369         (
370             "Foam::labelPair Foam::globalIndexAndTransform::encode"
371             "("
372                 "const label procI, "
373                 "const label index, "
374                 "const label transformIndex"
375             ")"
376         )
377             << "Overflow : encoding processor " << procI << " in base " << base_
378             << " exceeds capability of label (" << labelMax
379             << "). Please recompile with larger datatype for label."
380             << exit(FatalError);
381     }
383     return labelPair
384     (
385         index,
386         transformIndex + procI*base_
387     );
391 Foam::label Foam::globalIndexAndTransform::index
393     const labelPair& globalIAndTransform
396     return globalIAndTransform.first();
400 Foam::label Foam::globalIndexAndTransform::processor
402     const labelPair& globalIAndTransform
405     return globalIAndTransform.second()/base_;
409 Foam::label Foam::globalIndexAndTransform::transformIndex
411     const labelPair& globalIAndTransform
414     return globalIAndTransform.second() % base_;
418 Foam::label Foam::globalIndexAndTransform::nIndependentTransforms() const
420     return transforms_.size();
424 const Foam::List<Foam::vectorTensorTransform>&
425 Foam::globalIndexAndTransform::transforms() const
427     return transforms_;
431 const Foam::List<Foam::vectorTensorTransform>&
432 Foam::globalIndexAndTransform::transformPermutations() const
434     return transformPermutations_;
438 Foam::label Foam::globalIndexAndTransform::nullTransformIndex() const
440     return nullTransformIndex_;
444 const Foam::List<Foam::Pair<Foam::label> >&
445 Foam::globalIndexAndTransform::patchTransformSign() const
447     return patchTransformSign_;
451 const Foam::vectorTensorTransform& Foam::globalIndexAndTransform::transform
453     label transformIndex
454 ) const
456     return transformPermutations_[transformIndex];
460 Foam::labelList Foam::globalIndexAndTransform::transformIndicesForPatches
462     const labelHashSet& patchIs
463 ) const
465     List<label> permutation(transforms_.size(), 0);
467     labelList selectedTransformIs(0);
469     if (patchIs.empty() || transforms_.empty())
470     {
471         return selectedTransformIs;
472     }
474     forAllConstIter(labelHashSet, patchIs, iter)
475     {
476         label patchI = iter.key();
478         const Pair<label>& transSign = patchTransformSign_[patchI];
480         label matchTransI = transSign.first();
482         if (matchTransI > -1)
483         {
484             label sign = transSign.second();
486             // If this transform been found already by a patch?
487             if (permutation[matchTransI] != 0)
488             {
489                 // If so, if they have opposite signs, then this is
490                 // considered an error.  They are allowed to be the
491                 // same sign, but this only results in a single
492                 // transform.
493                 if (permutation[matchTransI] != sign)
494                 {
495                     FatalErrorIn
496                     (
497                         "const Foam::List<Foam::vectorTensorTransform>& "
498                         "Foam::globalIndexAndTransform::transformsForPatches"
499                         "("
500                             "const labelList& patchIs"
501                         ") const"
502                     )
503                         << "More than one patch accessing the same transform "
504                         << "but not of the same sign."
505                         << exit(FatalError);
506                 }
507             }
508             else
509             {
510                 permutation[matchTransI] = sign;
511             }
512         }
513     }
515     label nUsedTrans = round(sum(mag(permutation)));
517     if (nUsedTrans == 0)
518     {
519         return selectedTransformIs;
520     }
522     // Number of selected transformations
523     label nSelTrans = pow(2, nUsedTrans) - 1;
525     // Pout<< nl << permutation << nl << endl;
527     selectedTransformIs.setSize(nSelTrans);
529     switch (nUsedTrans)
530     {
531         case 1:
532         {
533             selectedTransformIs[0] = encodeTransformIndex(permutation);
535             break;
536         }
537         case 2:
538         {
539             List<label> tempPermutation = permutation;
541             label a = 0;
542             label b = 1;
544             // When there are two selected transforms out of three, we
545             // need to choose which of them are being permuted
546             if (transforms_.size() > nUsedTrans)
547             {
548                 if (permutation[0] == 0)
549                 {
550                     a = 1;
551                     b = 2;
552                 }
553                 else if (permutation[1] == 0)
554                 {
555                     a = 0;
556                     b = 2;
557                 }
558                 else if (permutation[2] == 0)
559                 {
560                     a = 0;
561                     b = 1;
562                 }
563             }
565             tempPermutation[a] = a;
566             tempPermutation[b] = permutation[b];
568             selectedTransformIs[0] = encodeTransformIndex(tempPermutation);
570             tempPermutation[a] = permutation[a];
571             tempPermutation[b] = a;
573             selectedTransformIs[1] = encodeTransformIndex(tempPermutation);
575             tempPermutation[a] = permutation[a];
576             tempPermutation[b] = permutation[b];
578             selectedTransformIs[2] = encodeTransformIndex(tempPermutation);
580             break;
581         }
582         case 3:
583         {
584             List<label> tempPermutation = permutation;
586             tempPermutation[0] = 0;
587             tempPermutation[1] = 0;
588             tempPermutation[2] = permutation[2];
590             selectedTransformIs[0] = encodeTransformIndex(tempPermutation);
592             tempPermutation[0] = 0;
593             tempPermutation[1] = permutation[1];
594             tempPermutation[2] = 0;
596             selectedTransformIs[1] = encodeTransformIndex(tempPermutation);
598             tempPermutation[0] = 0;
599             tempPermutation[1] = permutation[1];
600             tempPermutation[2] = permutation[2];
602             selectedTransformIs[2] = encodeTransformIndex(tempPermutation);
604             tempPermutation[0] = permutation[0];
605             tempPermutation[1] = 0;
606             tempPermutation[2] = 0;
608             selectedTransformIs[3] = encodeTransformIndex(tempPermutation);
610             tempPermutation[0] = permutation[0];
611             tempPermutation[1] = 0;
612             tempPermutation[2] = permutation[2];
614             selectedTransformIs[4] = encodeTransformIndex(tempPermutation);
616             tempPermutation[0] = permutation[0];
617             tempPermutation[1] = permutation[1];
618             tempPermutation[2] = 0;
620             selectedTransformIs[5] = encodeTransformIndex(tempPermutation);
622             tempPermutation[0] = permutation[0];
623             tempPermutation[1] = permutation[1];
624             tempPermutation[2] = permutation[2];
626             selectedTransformIs[6] = encodeTransformIndex(tempPermutation);
628             break;
629         }
630         default:
631         {
632             FatalErrorIn
633             (
634                 "const Foam::List<Foam::vectorTensorTransform>& "
635                 "Foam::globalIndexAndTransform::transformsForPatches"
636                 "("
637                     "const labelList& patchIs"
638                 ") const"
639             )
640                 << "Only 1-3 transforms are possible."
641                 << exit(FatalError);
642         }
643     }
645     return selectedTransformIs;
649 Foam::pointField Foam::globalIndexAndTransform::transformPatches
651     const labelHashSet& patchIs,
652     const point& pt
653 ) const
655     labelList transIs = transformIndicesForPatches(patchIs);
657     // Pout<< patchIs << nl << transIs << endl;
659     pointField transPts(transIs.size());
661     forAll(transIs, tII)
662     {
663         transPts[tII] = transformPermutations_[transIs[tII]].transformPosition
664         (
665             pt
666         );
667     }
669     return transPts;
673 // ************************************************************************* //