1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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 \*---------------------------------------------------------------------------*/
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 bool Foam::globalIndexAndTransform::less::operator()
36 label procA = globalIndexAndTransform::processor(a);
37 label procB = globalIndexAndTransform::processor(b);
43 else if (procA > procB)
50 label indexA = globalIndexAndTransform::index(a);
51 label indexB = globalIndexAndTransform::index(b);
57 else if (indexA > indexB)
64 label transformA = globalIndexAndTransform::transformIndex(a);
65 label transformB = globalIndexAndTransform::transformIndex(b);
67 return transformA < transformB;
73 Foam::label Foam::globalIndexAndTransform::encodeTransformIndex
75 const List<label>& permutationIndices
78 if (permutationIndices.size() != transforms_.size())
82 "Foam::label encodeTransformIndex"
84 "const List<label>& permutationIndices,"
87 << "permutationIndices " << permutationIndices
88 << "are of a different size to the number of independent transforms"
92 label transformIndex = 0;
96 forAll(transforms_, b)
98 if (mag(permutationIndices[b]) > 1)
102 "Foam::label encodeTransformIndex"
104 "const List<label>& permutationIndices,"
107 << "permutationIndices " << permutationIndices
108 << "are illegal, they must all be only -1, 0 or +1"
109 << abort(FatalError);
112 transformIndex += (permutationIndices[b] + 1)*w;
117 return transformIndex;
121 Foam::label Foam::globalIndexAndTransform::encodeTransformIndex
123 const FixedList<Foam::label, 3>& permutation
126 if (nIndependentTransforms() == 0)
130 if (nIndependentTransforms() == 1)
132 return permutation[0]+1;
134 else if (nIndependentTransforms() == 2)
136 return (permutation[1]+1)*3 + (permutation[0]+1);
142 + (permutation[1]+1)*3
143 + (permutation[0]+1);
148 Foam::FixedList<Foam::label, 3>
149 Foam::globalIndexAndTransform::decodeTransformIndex
151 const label transformIndex
154 FixedList<label, 3> permutation(0);
156 label t = transformIndex;
157 if (nIndependentTransforms() > 0)
159 permutation[0] = (t%3)-1;
160 if (nIndependentTransforms() > 1)
163 permutation[1] = (t%3)-1;
164 if (nIndependentTransforms() > 2)
167 permutation[2] = (t%3)-1;
178 "globalIndexAndTransform::decodeTransformIndex(const label)"
179 ) << "transformIndex : " << transformIndex
180 << " has more than 3 fields."
181 << abort(FatalError);
189 Foam::label Foam::globalIndexAndTransform::addToTransformIndex
191 const label transformIndex,
193 const bool isSendingSide
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)
204 FixedList<label, 3> permutation = decodeTransformIndex(transformIndex);
207 // Add patch transform
208 // ~~~~~~~~~~~~~~~~~~~
210 label sign = transSign.second();
217 // If this transform been found already by a patch?
218 if (permutation[matchTransI] != 0)
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
229 else if (sign == permutation[matchTransI])
234 "Foam::globalIndexAndTransform::addToTransformIndex\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
249 permutation[matchTransI] = 0;
254 permutation[matchTransI] = sign;
258 // Re-encode permutation
259 // ~~~~~~~~~~~~~~~~~~~~~
261 return encodeTransformIndex(permutation);
265 return transformIndex;
270 Foam::label Foam::globalIndexAndTransform::minimumTransformIndex
272 const label transformIndex0,
273 const label transformIndex1
276 if (transformIndex0 == transformIndex1)
278 return transformIndex0;
282 // Count number of transforms
283 FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
285 forAll(permutation0, i)
287 if (permutation0[i] != 0)
293 FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);
295 forAll(permutation1, i)
297 if (permutation1[i] != 0)
305 return transformIndex0;
309 return transformIndex1;
314 Foam::label Foam::globalIndexAndTransform::subtractTransformIndex
316 const label transformIndex0,
317 const label transformIndex1
320 FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
321 FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);
323 forAll(permutation0, i)
325 permutation0[i] -= permutation1[i];
328 return encodeTransformIndex(permutation0);
332 Foam::labelPair Foam::globalIndexAndTransform::encode
335 const label transformIndex
338 return encode(Pstream::myProcNo(), index, transformIndex);
342 Foam::labelPair Foam::globalIndexAndTransform::encode
346 const label transformIndex
349 if (transformIndex < 0 || transformIndex >= base_)
353 "Foam::labelPair Foam::globalIndexAndTransform::encode"
355 "const label procI, "
356 "const label index, "
357 "const label transformIndex"
360 << "TransformIndex " << transformIndex
361 << " is outside allowed range of 0 to "
363 << abort(FatalError);
366 if (procI > labelMax/base_)
370 "Foam::labelPair Foam::globalIndexAndTransform::encode"
372 "const label procI, "
373 "const label index, "
374 "const label transformIndex"
377 << "Overflow : encoding processor " << procI << " in base " << base_
378 << " exceeds capability of label (" << labelMax
379 << "). Please recompile with larger datatype for label."
386 transformIndex + procI*base_
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
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
456 return transformPermutations_[transformIndex];
460 Foam::labelList Foam::globalIndexAndTransform::transformIndicesForPatches
462 const labelHashSet& patchIs
465 List<label> permutation(transforms_.size(), 0);
467 labelList selectedTransformIs(0);
469 if (patchIs.empty() || transforms_.empty())
471 return selectedTransformIs;
474 forAllConstIter(labelHashSet, patchIs, iter)
476 label patchI = iter.key();
478 const Pair<label>& transSign = patchTransformSign_[patchI];
480 label matchTransI = transSign.first();
482 if (matchTransI > -1)
484 label sign = transSign.second();
486 // If this transform been found already by a patch?
487 if (permutation[matchTransI] != 0)
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
493 if (permutation[matchTransI] != sign)
497 "const Foam::List<Foam::vectorTensorTransform>& "
498 "Foam::globalIndexAndTransform::transformsForPatches"
500 "const labelList& patchIs"
503 << "More than one patch accessing the same transform "
504 << "but not of the same sign."
510 permutation[matchTransI] = sign;
515 label nUsedTrans = round(sum(mag(permutation)));
519 return selectedTransformIs;
522 // Number of selected transformations
523 label nSelTrans = pow(2, nUsedTrans) - 1;
525 // Pout<< nl << permutation << nl << endl;
527 selectedTransformIs.setSize(nSelTrans);
533 selectedTransformIs[0] = encodeTransformIndex(permutation);
539 List<label> tempPermutation = permutation;
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)
548 if (permutation[0] == 0)
553 else if (permutation[1] == 0)
558 else if (permutation[2] == 0)
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);
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);
634 "const Foam::List<Foam::vectorTensorTransform>& "
635 "Foam::globalIndexAndTransform::transformsForPatches"
637 "const labelList& patchIs"
640 << "Only 1-3 transforms are possible."
645 return selectedTransformIs;
649 Foam::pointField Foam::globalIndexAndTransform::transformPatches
651 const labelHashSet& patchIs,
655 labelList transIs = transformIndicesForPatches(patchIs);
657 // Pout<< patchIs << nl << transIs << endl;
659 pointField transPts(transIs.size());
663 transPts[tII] = transformPermutations_[transIs[tII]].transformPosition
673 // ************************************************************************* //