globalIndexAndTransformI.H 16.7 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2010-2011 OpenCFD Ltd.
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

26
#include "polyMesh.H"
27

28
// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
29

30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
bool Foam::globalIndexAndTransform::less::operator()
(
    const labelPair& a,
    const labelPair& b
) const
{
    label procA = globalIndexAndTransform::processor(a);
    label procB = globalIndexAndTransform::processor(b);

    if (procA < procB)
    {
        return true;
    }
    else if (procA > procB)
    {
        return false;
    }
    else
    {
        // Equal proc.
        label indexA = globalIndexAndTransform::index(a);
        label indexB = globalIndexAndTransform::index(b);

        if (indexA < indexB)
        {
            return true;
        }
        else if (indexA > indexB)
        {
            return false;
        }
        else
        {
            // Equal index
            label transformA = globalIndexAndTransform::transformIndex(a);
            label transformB = globalIndexAndTransform::transformIndex(b);

            return transformA < transformB;
        }
    }
}


73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
Foam::label Foam::globalIndexAndTransform::encodeTransformIndex
(
    const List<label>& permutationIndices
) const
{
    if (permutationIndices.size() != transforms_.size())
    {
        FatalErrorIn
        (
            "Foam::label encodeTransformIndex"
            "("
                "const List<label>& permutationIndices,"
            ") const"
        )
            << "permutationIndices " << permutationIndices
            << "are of a different size to the number of independent transforms"
            << abort(FatalError);
    }

    label transformIndex = 0;

    label w = 1;

96
    forAll(transforms_, b)
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
    {
        if (mag(permutationIndices[b]) > 1)
        {
            FatalErrorIn
            (
                "Foam::label encodeTransformIndex"
                "("
                "const List<label>& permutationIndices,"
                ") const"
            )
                << "permutationIndices " << permutationIndices
                << "are illegal, they must all be only -1, 0 or +1"
                << abort(FatalError);
        }

        transformIndex += (permutationIndices[b] + 1)*w;

        w *= 3;
    }

    return transformIndex;
}


121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
Foam::label Foam::globalIndexAndTransform::encodeTransformIndex
(
    const FixedList<Foam::label, 3>& permutation
) const
{
    if (nIndependentTransforms() == 0)
    {
        return 0;
    }
    if (nIndependentTransforms() == 1)
    {
        return permutation[0]+1;
    }
    else if (nIndependentTransforms() == 2)
    {
        return (permutation[1]+1)*3 + (permutation[0]+1);
    }
    else
    {
        return
            (permutation[2]+1)*9
          + (permutation[1]+1)*3
          + (permutation[0]+1);
    }
}


148
149
150
151
Foam::FixedList<Foam::label, 3>
Foam::globalIndexAndTransform::decodeTransformIndex
(
    const label transformIndex
152
) const
153
{
154
    FixedList<label, 3> permutation(0);
155
156

    label t = transformIndex;
157
158
159
160
161
162
163
164
165
166
167
168
169
170
    if (nIndependentTransforms() > 0)
    {
        permutation[0] = (t%3)-1;
        if (nIndependentTransforms() > 1)
        {
            t /= 3;
            permutation[1] = (t%3)-1;
            if (nIndependentTransforms() > 2)
            {
                t /= 3;
                permutation[2] = (t%3)-1;
            }
        }
    }
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188

#   ifdef FULLDEBUG
    t /= 3;
    if (t != 0)
    {
        FatalErrorIn
        (
            "globalIndexAndTransform::decodeTransformIndex(const label)"
        )   << "transformIndex : " << transformIndex
            << " has more than 3 fields."
            << abort(FatalError);
    }
#   endif

    return permutation;
}


189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
Foam::label Foam::globalIndexAndTransform::addToTransformIndex
(
    const label transformIndex,
    const label patchI,
    const bool isSendingSide
) const
{
    const Pair<label>& transSign = patchTransformSign_[patchI];

    label matchTransI = transSign.first();

    // Hardcoded for max 3 transforms only!

    if (matchTransI > -1 && matchTransI < 3)
    {
204
        FixedList<label, 3> permutation = decodeTransformIndex(transformIndex);
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219


        // Add patch transform
        // ~~~~~~~~~~~~~~~~~~~

        label sign = transSign.second();
        if (!isSendingSide)
        {
            sign = -sign;
        }


        // If this transform been found already by a patch?
        if (permutation[matchTransI] != 0)
        {
220
221
222
223
224
225
226
227
228
229
            if (sign == 0)
            {
                // sent from patch without a transformation. Do nothing.
                FatalErrorIn("globalIndexAndTransform::addToTransformIndex(..)")
                    << "patch:" << mesh_.boundaryMesh()[patchI].name()
                    << " transform:" << matchTransI << " sign:" << sign
                    << "  current transforms:" << permutation
                    << exit(FatalError);
            }
            else if (sign == permutation[matchTransI])
230
231
232
233
234
235
            {
                FatalErrorIn
                (
                    "Foam::label "
                    "Foam::globalIndexAndTransform::addToTransformIndex\n"
                    "(\n"
236
237
238
                        "const label,\n"
                        "const label,\n"
                        "const bool\n"
239
240
241
242
243
                    ") const\n"
                )   << "More than one patch accessing the same transform "
                    << "but not of the same sign." << endl
                    << "patch:" << mesh_.boundaryMesh()[patchI].name()
                    << " transform:" << matchTransI << " sign:" << sign
244
                    << "  current transforms:" << permutation
245
246
                    << exit(FatalError);
            }
247
248
249
250
            else
            {
                permutation[matchTransI] = 0;
            }
251
252
253
254
255
256
257
        }
        else
        {
            permutation[matchTransI] = sign;
        }


258
259
        // Re-encode permutation
        // ~~~~~~~~~~~~~~~~~~~~~
260

261
        return encodeTransformIndex(permutation);
262
263
264
265
266
267
268
269
    }
    else
    {
        return transformIndex;
    }
}


270
Foam::label Foam::globalIndexAndTransform::minimumTransformIndex
271
272
273
(
    const label transformIndex0,
    const label transformIndex1
274
) const
275
{
276
    if (transformIndex0 == transformIndex1)
277
    {
278
        return transformIndex0;
279
280
281
    }


282
    // Count number of transforms
283
    FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
284
    label n0 = 0;
285
286
    forAll(permutation0, i)
    {
287
        if (permutation0[i] != 0)
288
        {
289
            n0++;
290
        }
291
292
293
294
295
296
297
    }

    FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);
    label n1 = 0;
    forAll(permutation1, i)
    {
        if (permutation1[i] != 0)
298
        {
299
            n1++;
300
301
        }
    }
302
303
304
305
306
307
308
309
310

    if (n0 <= n1)
    {
        return transformIndex0;
    }
    else
    {
        return transformIndex1;
    }
311
312
313
}


314
315
316
317
Foam::label Foam::globalIndexAndTransform::subtractTransformIndex
(
    const label transformIndex0,
    const label transformIndex1
318
) const
319
320
321
322
323
324
325
326
327
{
    FixedList<label, 3> permutation0 = decodeTransformIndex(transformIndex0);
    FixedList<label, 3> permutation1 = decodeTransformIndex(transformIndex1);

    forAll(permutation0, i)
    {
        permutation0[i] -= permutation1[i];
    }

328
    return encodeTransformIndex(permutation0);
329
330
331
}


332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
Foam::labelPair Foam::globalIndexAndTransform::encode
(
    const label index,
    const label transformIndex
)
{
    return encode(Pstream::myProcNo(), index, transformIndex);
}


Foam::labelPair Foam::globalIndexAndTransform::encode
(
    const label procI,
    const label index,
    const label transformIndex
)
{
    if (transformIndex < 0 || transformIndex >= base_)
    {
        FatalErrorIn
        (
            "Foam::labelPair Foam::globalIndexAndTransform::encode"
            "("
                "const label procI, "
                "const label index, "
                "const label transformIndex"
            ")"
        )
            << "TransformIndex " << transformIndex
            << " is outside allowed range of 0 to "
            << base_ - 1
            << abort(FatalError);
    }

    if (procI > labelMax/base_)
    {
        FatalErrorIn
        (
            "Foam::labelPair Foam::globalIndexAndTransform::encode"
            "("
                "const label procI, "
                "const label index, "
                "const label transformIndex"
            ")"
        )
            << "Overflow : encoding processor " << procI << " in base " << base_
            << " exceeds capability of label (" << labelMax
            << "). Please recompile with larger datatype for label."
            << exit(FatalError);
    }

    return labelPair
    (
        index,
        transformIndex + procI*base_
    );
}


391
Foam::label Foam::globalIndexAndTransform::index
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
(
    const labelPair& globalIAndTransform
)
{
    return globalIAndTransform.first();
}


Foam::label Foam::globalIndexAndTransform::processor
(
    const labelPair& globalIAndTransform
)
{
    return globalIAndTransform.second()/base_;
}


Foam::label Foam::globalIndexAndTransform::transformIndex
(
    const labelPair& globalIAndTransform
)
{
    return globalIAndTransform.second() % base_;
}


Foam::label Foam::globalIndexAndTransform::nIndependentTransforms() const
{
    return transforms_.size();
}


424
const Foam::List<Foam::vectorTensorTransform>&
425
426
427
428
429
430
Foam::globalIndexAndTransform::transforms() const
{
    return transforms_;
}


431
const Foam::List<Foam::vectorTensorTransform>&
432
433
434
435
436
437
Foam::globalIndexAndTransform::transformPermutations() const
{
    return transformPermutations_;
}


438
439
440
441
442
443
Foam::label Foam::globalIndexAndTransform::nullTransformIndex() const
{
    return nullTransformIndex_;
}


444
445
446
447
448
449
450
const Foam::List<Foam::Pair<Foam::label> >&
Foam::globalIndexAndTransform::patchTransformSign() const
{
    return patchTransformSign_;
}


451
const Foam::vectorTensorTransform& Foam::globalIndexAndTransform::transform
452
453
454
455
456
457
458
459
(
    label transformIndex
) const
{
    return transformPermutations_[transformIndex];
}


460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
Foam::labelList Foam::globalIndexAndTransform::transformIndicesForPatches
(
    const labelHashSet& patchIs
) const
{
    List<label> permutation(transforms_.size(), 0);

    labelList selectedTransformIs(0);

    if (patchIs.empty() || transforms_.empty())
    {
        return selectedTransformIs;
    }

    forAllConstIter(labelHashSet, patchIs, iter)
    {
        label patchI = iter.key();

        const Pair<label>& transSign = patchTransformSign_[patchI];

        label matchTransI = transSign.first();

        if (matchTransI > -1)
        {
            label sign = transSign.second();

            // If this transform been found already by a patch?
            if (permutation[matchTransI] != 0)
            {
                // If so, if they have opposite signs, then this is
                // considered an error.  They are allowed to be the
                // same sign, but this only results in a single
                // transform.
                if (permutation[matchTransI] != sign)
                {
                    FatalErrorIn
                    (
                        "const Foam::List<Foam::vectorTensorTransform>& "
                        "Foam::globalIndexAndTransform::transformsForPatches"
                        "("
                            "const labelList& patchIs"
                        ") const"
                    )
                        << "More than one patch accessing the same transform "
                        << "but not of the same sign."
                        << exit(FatalError);
                }
            }
            else
            {
                permutation[matchTransI] = sign;
            }
        }
    }

mattijs's avatar
mattijs committed
515
    label nUsedTrans = round(sum(mag(permutation)));
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662

    if (nUsedTrans == 0)
    {
        return selectedTransformIs;
    }

    // Number of selected transformations
    label nSelTrans = pow(2, nUsedTrans) - 1;

    // Pout<< nl << permutation << nl << endl;

    selectedTransformIs.setSize(nSelTrans);

    switch (nUsedTrans)
    {
        case 1:
        {
            selectedTransformIs[0] = encodeTransformIndex(permutation);

            break;
        }
        case 2:
        {
            List<label> tempPermutation = permutation;

            label a = 0;
            label b = 1;

            // When there are two selected transforms out of three, we
            // need to choose which of them are being permuted
            if (transforms_.size() > nUsedTrans)
            {
                if (permutation[0] == 0)
                {
                    a = 1;
                    b = 2;
                }
                else if (permutation[1] == 0)
                {
                    a = 0;
                    b = 2;
                }
                else if (permutation[2] == 0)
                {
                    a = 0;
                    b = 1;
                }
            }

            tempPermutation[a] = a;
            tempPermutation[b] = permutation[b];

            selectedTransformIs[0] = encodeTransformIndex(tempPermutation);

            tempPermutation[a] = permutation[a];
            tempPermutation[b] = a;

            selectedTransformIs[1] = encodeTransformIndex(tempPermutation);

            tempPermutation[a] = permutation[a];
            tempPermutation[b] = permutation[b];

            selectedTransformIs[2] = encodeTransformIndex(tempPermutation);

            break;
        }
        case 3:
        {
            List<label> tempPermutation = permutation;

            tempPermutation[0] = 0;
            tempPermutation[1] = 0;
            tempPermutation[2] = permutation[2];

            selectedTransformIs[0] = encodeTransformIndex(tempPermutation);

            tempPermutation[0] = 0;
            tempPermutation[1] = permutation[1];
            tempPermutation[2] = 0;

            selectedTransformIs[1] = encodeTransformIndex(tempPermutation);

            tempPermutation[0] = 0;
            tempPermutation[1] = permutation[1];
            tempPermutation[2] = permutation[2];

            selectedTransformIs[2] = encodeTransformIndex(tempPermutation);

            tempPermutation[0] = permutation[0];
            tempPermutation[1] = 0;
            tempPermutation[2] = 0;

            selectedTransformIs[3] = encodeTransformIndex(tempPermutation);

            tempPermutation[0] = permutation[0];
            tempPermutation[1] = 0;
            tempPermutation[2] = permutation[2];

            selectedTransformIs[4] = encodeTransformIndex(tempPermutation);

            tempPermutation[0] = permutation[0];
            tempPermutation[1] = permutation[1];
            tempPermutation[2] = 0;

            selectedTransformIs[5] = encodeTransformIndex(tempPermutation);

            tempPermutation[0] = permutation[0];
            tempPermutation[1] = permutation[1];
            tempPermutation[2] = permutation[2];

            selectedTransformIs[6] = encodeTransformIndex(tempPermutation);

            break;
        }
        default:
        {
            FatalErrorIn
            (
                "const Foam::List<Foam::vectorTensorTransform>& "
                "Foam::globalIndexAndTransform::transformsForPatches"
                "("
                    "const labelList& patchIs"
                ") const"
            )
                << "Only 1-3 transforms are possible."
                << exit(FatalError);
        }
    }

    return selectedTransformIs;
}


Foam::pointField Foam::globalIndexAndTransform::transformPatches
(
    const labelHashSet& patchIs,
    const point& pt
) const
{
    labelList transIs = transformIndicesForPatches(patchIs);

    // Pout<< patchIs << nl << transIs << endl;

    pointField transPts(transIs.size());

    forAll(transIs, tII)
    {
663
664
665
666
        transPts[tII] = transformPermutations_[transIs[tII]].transformPosition
        (
            pt
        );
667
668
669
670
671
672
    }

    return transPts;
}


673
// ************************************************************************* //