AMIInterpolation.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "AMIInterpolation.H"
30#include "meshTools.H"
31#include "mapDistribute.H"
32#include "flipOp.H"
33#include "profiling.H"
34#include "triPointRef.H"
35#include "OFstream.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
44}
45
47
48
49// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50
53(
54 const primitivePatch& patch
55) const
56{
57 treeBoundBox bb(patch.points(), patch.meshPoints());
58 bb.inflate(0.01);
59
61 (
63 (
64 false,
65 patch,
67 ),
68 bb, // overall search domain
69 8, // maxLevel
70 10, // leaf size
71 3.0 // duplicity
72 );
73}
74
75
77(
78 const primitivePatch& srcPatch,
79 const primitivePatch& tgtPatch
80) const
81{
82 // Either not parallel or no faces on any processor
83 label proci = 0;
84
85 if (Pstream::parRun())
86 {
87 const bitSet hasFaces
88 (
89 UPstream::listGatherValues<bool>
90 (
91 srcPatch.size() > 0 || tgtPatch.size() > 0
92 )
93 );
94
95 const auto nHaveFaces = hasFaces.count();
96
97 if (nHaveFaces == 1)
98 {
99 proci = hasFaces.find_first();
101 << "AMI local to processor" << proci << endl;
102 }
103 else if (nHaveFaces > 1)
104 {
105 proci = -1;
107 << "AMI split across multiple processors" << endl;
108 }
109
110 Pstream::broadcast(proci);
111 }
112
113 return proci;
114}
115
116
118(
119 const searchableSurface& surf,
120 pointField& pts
121) const
122{
123 addProfiling(ami, "AMIInterpolation::projectPointsToSurface");
124
125 DebugInfo<< "AMI: projecting points to surface" << endl;
126
127 List<pointIndexHit> nearInfo;
128
129 surf.findNearest(pts, scalarField(pts.size(), GREAT), nearInfo);
130
131 label nMiss = 0;
132 forAll(nearInfo, i)
133 {
134 const pointIndexHit& pi = nearInfo[i];
135
136 if (pi.hit())
137 {
138 pts[i] = pi.hitPoint();
139 }
140 else
141 {
142 // Point remains unchanged
143 ++nMiss;
144 }
145 }
146
147 if (nMiss > 0)
148 {
150 << "Error projecting points to surface: "
151 << nMiss << " faces could not be determined"
152 << abort(FatalError);
153 }
154}
155
156
158(
159 const scalarList& patchAreas,
160 const word& patchName,
161 const labelListList& addr,
162 scalarListList& wght,
163 scalarField& wghtSum,
164 const bool conformal,
165 const bool output,
166 const scalar lowWeightTol
167)
168{
169 addProfiling(ami, "AMIInterpolation::normaliseWeights");
170
171 // Normalise the weights
172 wghtSum.setSize(wght.size(), 0.0);
173 label nLowWeight = 0;
174
175 forAll(wght, facei)
176 {
177 scalarList& w = wght[facei];
178
179 if (w.size())
180 {
181 scalar denom = patchAreas[facei];
182
183 scalar s = sum(w);
184 scalar t = s/denom;
185 if (conformal)
186 {
187 denom = s;
188 }
189
190 forAll(w, i)
191 {
192 w[i] /= denom;
193 }
194
195 wghtSum[facei] = t;
196 if (t < lowWeightTol)
197 {
198 ++nLowWeight;
199 }
200 }
201 else
202 {
203 wghtSum[facei] = 0;
204 }
205 }
206
207 if (output)
208 {
209 const label nFace = returnReduce(wght.size(), sumOp<label>());
210
211 if (nFace)
212 {
213 Info<< indent
214 << "AMI: Patch " << patchName
215 << " sum(weights)"
216 << " min:" << gMin(wghtSum)
217 << " max:" << gMax(wghtSum)
218 << " average:" << gAverage(wghtSum) << nl;
219
220 const label nLow = returnReduce(nLowWeight, sumOp<label>());
221
222 if (nLow)
223 {
224 Info<< indent
225 << "AMI: Patch " << patchName
226 << " identified " << nLow
227 << " faces with weights less than " << lowWeightTol
228 << endl;
229 }
230 }
231 }
232}
233
234
236(
237 const autoPtr<mapDistribute>& targetMapPtr,
238 const scalarList& fineSrcMagSf,
239 const labelListList& fineSrcAddress,
240 const scalarListList& fineSrcWeights,
241
242 const labelList& sourceRestrictAddressing,
243 const labelList& targetRestrictAddressing,
244
245 scalarList& srcMagSf,
246 labelListList& srcAddress,
247 scalarListList& srcWeights,
248 scalarField& srcWeightsSum,
250)
251{
252 addProfiling(ami, "AMIInterpolation::agglomerate");
253
254 label sourceCoarseSize =
255 (
256 sourceRestrictAddressing.size()
257 ? max(sourceRestrictAddressing)+1
258 : 0
259 );
260
261 label targetCoarseSize =
262 (
263 targetRestrictAddressing.size()
264 ? max(targetRestrictAddressing)+1
265 : 0
266 );
267
268 // Agglomerate face areas
269 {
270 srcMagSf.setSize(sourceRestrictAddressing.size(), 0.0);
271 forAll(sourceRestrictAddressing, facei)
272 {
273 label coarseFacei = sourceRestrictAddressing[facei];
274 srcMagSf[coarseFacei] += fineSrcMagSf[facei];
275 }
276 }
277
278 // Agglomerate weights and indices
279 if (targetMapPtr)
280 {
281 const mapDistribute& map = *targetMapPtr;
282
283 // Get all restriction addressing.
284 labelList allRestrict(targetRestrictAddressing);
285 map.distribute(allRestrict);
286
287 // So now we have agglomeration of the target side in
288 // allRestrict:
289 // 0..size-1 : local agglomeration (= targetRestrictAddressing
290 // (but potentially permutated))
291 // size.. : agglomeration data from other processors
292
293
294 // The trickiness in this algorithm is finding out the compaction
295 // of the remote data (i.e. allocation of the coarse 'slots'). We could
296 // either send across the slot compaction maps or just make sure
297 // that we allocate the slots in exactly the same order on both sending
298 // and receiving side (e.g. if the submap is set up to send 4 items,
299 // the constructMap is also set up to receive 4 items.
300
301
302 // Short note about the various types of indices:
303 // - face indices : indices into the geometry.
304 // - coarse face indices : how the faces get agglomerated
305 // - transferred data : how mapDistribute sends/receives data,
306 // - slots : indices into data after distribution (e.g. stencil,
307 // srcAddress/tgtAddress). Note: for fully local addressing
308 // the slots are equal to face indices.
309 // A mapDistribute has:
310 // - a subMap : these are face indices
311 // - a constructMap : these are from 'transferred-data' to slots
312
313 labelListList tgtSubMap(Pstream::nProcs());
314
315 // Local subMap is just identity
316 {
317 tgtSubMap[Pstream::myProcNo()] = identity(targetCoarseSize);
318 }
319
320 forAll(map.subMap(), proci)
321 {
322 if (proci != Pstream::myProcNo())
323 {
324 // Combine entries that point to the same coarse element.
325 // The important bit is to loop over the data (and hand out
326 // compact indices ) in 'transferred data' order. This
327 // guarantees that we're doing exactly the
328 // same on sending and receiving side - e.g. the fourth element
329 // in the subMap is the fourth element received in the
330 // constructMap
331
332 const labelList& elems = map.subMap()[proci];
333 const labelList& elemsMap =
335 labelList& newSubMap = tgtSubMap[proci];
336 newSubMap.setSize(elems.size());
337
338 labelList oldToNew(targetCoarseSize, -1);
339 label newi = 0;
340
341 for (const label elemi : elems)
342 {
343 label fineElem = elemsMap[elemi];
344 label coarseElem = allRestrict[fineElem];
345 if (oldToNew[coarseElem] == -1)
346 {
347 oldToNew[coarseElem] = newi;
348 newSubMap[newi] = coarseElem;
349 ++newi;
350 }
351 }
352 newSubMap.setSize(newi);
353 }
354 }
355
356 // Reconstruct constructMap by combining entries. Note that order
357 // of handing out indices should be the same as loop above to compact
358 // the sending map
359
360 labelListList tgtConstructMap(Pstream::nProcs());
361
362 // Local constructMap is just identity
363 {
364 tgtConstructMap[Pstream::myProcNo()] = identity(targetCoarseSize);
365 }
366
367 labelList tgtCompactMap(map.constructSize());
368
369 {
370 // Note that in special cases (e.g. 'appending' two AMIs) the
371 // local size after distributing can be longer than the number
372 // of faces. I.e. it duplicates elements.
373 // Since we don't know this size instead we loop over all
374 // reachable elements (using the local constructMap)
375
376 const labelList& elemsMap = map.constructMap()[Pstream::myProcNo()];
377 for (const label fineElem : elemsMap)
378 {
379 label coarseElem = allRestrict[fineElem];
380 tgtCompactMap[fineElem] = coarseElem;
381 }
382 }
383
384 label compacti = targetCoarseSize;
385
386 // Compact data from other processors
387 forAll(map.constructMap(), proci)
388 {
389 if (proci != Pstream::myProcNo())
390 {
391 // Combine entries that point to the same coarse element. All
392 // elements now are remote data so we cannot use any local
393 // data here - use allRestrict instead.
394 const labelList& elems = map.constructMap()[proci];
395
396 labelList& newConstructMap = tgtConstructMap[proci];
397 newConstructMap.setSize(elems.size());
398
399 if (elems.size())
400 {
401 // Get the maximum target coarse size for this set of
402 // received data.
403 label remoteTargetCoarseSize = labelMin;
404 for (const label elemi : elems)
405 {
406 remoteTargetCoarseSize = max
407 (
408 remoteTargetCoarseSize,
409 allRestrict[elemi]
410 );
411 }
412 remoteTargetCoarseSize += 1;
413
414 // Combine locally data coming from proci
415 labelList oldToNew(remoteTargetCoarseSize, -1);
416 label newi = 0;
417
418 for (const label fineElem : elems)
419 {
420 // fineElem now points to section from proci
421 label coarseElem = allRestrict[fineElem];
422 if (oldToNew[coarseElem] == -1)
423 {
424 oldToNew[coarseElem] = newi;
425 tgtCompactMap[fineElem] = compacti;
426 newConstructMap[newi] = compacti++;
427 ++newi;
428 }
429 else
430 {
431 // Get compact index
432 label compacti = oldToNew[coarseElem];
433 tgtCompactMap[fineElem] = newConstructMap[compacti];
434 }
435 }
436 newConstructMap.setSize(newi);
437 }
438 }
439 }
440
441 srcAddress.setSize(sourceCoarseSize);
442 srcWeights.setSize(sourceCoarseSize);
443
444 forAll(fineSrcAddress, facei)
445 {
446 // All the elements contributing to facei. Are slots in
447 // mapDistribute'd data.
448 const labelList& elems = fineSrcAddress[facei];
449 const scalarList& weights = fineSrcWeights[facei];
450 const scalar fineArea = fineSrcMagSf[facei];
451
452 label coarseFacei = sourceRestrictAddressing[facei];
453
454 labelList& newElems = srcAddress[coarseFacei];
455 scalarList& newWeights = srcWeights[coarseFacei];
456
457 forAll(elems, i)
458 {
459 label elemi = elems[i];
460 label coarseElemi = tgtCompactMap[elemi];
461
462 label index = newElems.find(coarseElemi);
463 if (index == -1)
464 {
465 newElems.append(coarseElemi);
466 newWeights.append(fineArea*weights[i]);
467 }
468 else
469 {
470 newWeights[index] += fineArea*weights[i];
471 }
472 }
473 }
474
475 tgtMap.reset
476 (
477 new mapDistribute
478 (
479 compacti,
480 std::move(tgtSubMap),
481 std::move(tgtConstructMap)
482 )
483 );
484 }
485 else
486 {
487 srcAddress.setSize(sourceCoarseSize);
488 srcWeights.setSize(sourceCoarseSize);
489
490 forAll(fineSrcAddress, facei)
491 {
492 // All the elements contributing to facei. Are slots in
493 // mapDistribute'd data.
494 const labelList& elems = fineSrcAddress[facei];
495 const scalarList& weights = fineSrcWeights[facei];
496 const scalar fineArea = fineSrcMagSf[facei];
497
498 label coarseFacei = sourceRestrictAddressing[facei];
499
500 labelList& newElems = srcAddress[coarseFacei];
501 scalarList& newWeights = srcWeights[coarseFacei];
502
503 forAll(elems, i)
504 {
505 const label elemi = elems[i];
506 const label coarseElemi = targetRestrictAddressing[elemi];
507
508 const label index = newElems.find(coarseElemi);
509 if (index == -1)
510 {
511 newElems.append(coarseElemi);
512 newWeights.append(fineArea*weights[i]);
513 }
514 else
515 {
516 newWeights[index] += fineArea*weights[i];
517 }
518 }
519 }
520 }
521
522 // Weights normalisation
523 normaliseWeights
524 (
525 srcMagSf,
526 "source",
527 srcAddress,
528 srcWeights,
529 srcWeightsSum,
530 true,
531 false,
532 -1
533 );
534}
535
536
537// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
538
540(
541 const dictionary& dict,
542 const bool reverseTarget
543)
544:
545 requireMatch_(dict.getOrDefault("requireMatch", true)),
546 reverseTarget_(dict.getOrDefault("reverseTarget", reverseTarget)),
547 lowWeightCorrection_(dict.getOrDefault<scalar>("lowWeightCorrection", -1)),
548 singlePatchProc_(-999),
549 srcMagSf_(),
550 srcAddress_(),
551 srcWeights_(),
552 srcWeightsSum_(),
553 srcCentroids_(),
554 srcMapPtr_(nullptr),
555 tgtMagSf_(),
556 tgtAddress_(),
557 tgtWeights_(),
558 tgtWeightsSum_(),
559 tgtCentroids_(),
560 tgtMapPtr_(nullptr),
561 upToDate_(false)
562{}
563
564
566(
567 const bool requireMatch,
568 const bool reverseTarget,
569 const scalar lowWeightCorrection
570)
571:
572 requireMatch_(requireMatch),
573 reverseTarget_(reverseTarget),
574 lowWeightCorrection_(lowWeightCorrection),
575 singlePatchProc_(-999),
576 srcMagSf_(),
577 srcAddress_(),
578 srcWeights_(),
579 srcWeightsSum_(),
580 srcCentroids_(),
581 srcPatchPts_(),
582 srcMapPtr_(nullptr),
583 tgtMagSf_(),
584 tgtAddress_(),
585 tgtWeights_(),
586 tgtWeightsSum_(),
587 tgtCentroids_(),
588 tgtPatchPts_(),
589 tgtMapPtr_(nullptr),
590 upToDate_(false)
591{}
592
593
595(
596 const AMIInterpolation& fineAMI,
597 const labelList& sourceRestrictAddressing,
598 const labelList& targetRestrictAddressing
599)
600:
601 requireMatch_(fineAMI.requireMatch_),
602 reverseTarget_(fineAMI.reverseTarget_),
603 lowWeightCorrection_(-1.0),
604 singlePatchProc_(fineAMI.singlePatchProc_),
605 srcMagSf_(),
606 srcAddress_(),
607 srcWeights_(),
608 srcWeightsSum_(),
609 srcPatchPts_(),
610 srcMapPtr_(nullptr),
611 tgtMagSf_(),
612 tgtAddress_(),
613 tgtWeights_(),
614 tgtWeightsSum_(),
615 tgtPatchPts_(),
616 tgtMapPtr_(nullptr),
617 upToDate_(false)
618{
619 label sourceCoarseSize =
620 (
621 sourceRestrictAddressing.size()
622 ? max(sourceRestrictAddressing)+1
623 : 0
624 );
625
626 label neighbourCoarseSize =
627 (
628 targetRestrictAddressing.size()
629 ? max(targetRestrictAddressing)+1
630 : 0
631 );
632
633 if (debug & 2)
634 {
635 Pout<< "AMI: Creating addressing and weights as agglomeration of AMI :"
636 << " source:" << fineAMI.srcAddress().size()
637 << " target:" << fineAMI.tgtAddress().size()
638 << " coarse source size:" << sourceCoarseSize
639 << " neighbour source size:" << neighbourCoarseSize
640 << endl;
641 }
642
643 if
644 (
645 fineAMI.srcAddress().size() != sourceRestrictAddressing.size()
646 || fineAMI.tgtAddress().size() != targetRestrictAddressing.size()
647 )
648 {
650 << "Size mismatch." << nl
651 << "Source patch size:" << fineAMI.srcAddress().size() << nl
652 << "Source agglomeration size:"
653 << sourceRestrictAddressing.size() << nl
654 << "Target patch size:" << fineAMI.tgtAddress().size() << nl
655 << "Target agglomeration size:"
656 << targetRestrictAddressing.size()
657 << exit(FatalError);
658 }
659
660
661 // Agglomerate addresses and weights
662
664 (
665 fineAMI.tgtMapPtr_,
666 fineAMI.srcMagSf(),
667 fineAMI.srcAddress(),
668 fineAMI.srcWeights(),
669
670 sourceRestrictAddressing,
671 targetRestrictAddressing,
672
673 srcMagSf_,
678 );
679
681 (
682 fineAMI.srcMapPtr_,
683 fineAMI.tgtMagSf(),
684 fineAMI.tgtAddress(),
685 fineAMI.tgtWeights(),
686
687 targetRestrictAddressing,
688 sourceRestrictAddressing,
689
690 tgtMagSf_,
695 );
696}
697
698
700:
701 requireMatch_(ami.requireMatch_),
702 reverseTarget_(ami.reverseTarget_),
703 lowWeightCorrection_(ami.lowWeightCorrection_),
704 singlePatchProc_(ami.singlePatchProc_),
705 srcMagSf_(ami.srcMagSf_),
706 srcAddress_(ami.srcAddress_),
707 srcWeights_(ami.srcWeights_),
708 srcWeightsSum_(ami.srcWeightsSum_),
709 srcCentroids_(ami.srcCentroids_),
710 srcMapPtr_(nullptr),
711 tgtMagSf_(ami.tgtMagSf_),
712 tgtAddress_(ami.tgtAddress_),
713 tgtWeights_(ami.tgtWeights_),
714 tgtWeightsSum_(ami.tgtWeightsSum_),
715 tgtCentroids_(ami.tgtCentroids_),
716 tgtMapPtr_(nullptr),
717 upToDate_(false)
718{}
719
720
721// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
722
724(
725 const primitivePatch& srcPatch,
726 const primitivePatch& tgtPatch,
727 const autoPtr<searchableSurface>& surfPtr
728)
729{
730 if (upToDate_)
731 {
732 return false;
733 }
734
735 addProfiling(ami, "AMIInterpolation::calculate");
736
737 if (surfPtr)
738 {
739 srcPatchPts_ = srcPatch.points();
740 projectPointsToSurface(surfPtr(), srcPatchPts_);
741 tsrcPatch0_ = refPtr<primitivePatch>::New
742 (
743 SubList<face>(srcPatch),
744 srcPatchPts_
745 );
746
747 tgtPatchPts_ = tgtPatch.points();
748 projectPointsToSurface(surfPtr(), tgtPatchPts_);
749 ttgtPatch0_ = refPtr<primitivePatch>::New
750 (
751 SubList<face>(tgtPatch),
752 tgtPatchPts_
753 );
754 }
755 else
756 {
757 tsrcPatch0_.cref(srcPatch);
758 ttgtPatch0_.cref(tgtPatch);
759 }
760
761 label srcTotalSize = returnReduce(srcPatch.size(), sumOp<label>());
762 label tgtTotalSize = returnReduce(tgtPatch.size(), sumOp<label>());
763
764 if (srcTotalSize == 0)
765 {
766 DebugInfo<< "AMI: no source faces present - no addressing constructed"
767 << endl;
768
769 return false;
770 }
771
772 Info<< indent
773 << "AMI: Creating addressing and weights between "
774 << srcTotalSize << " source faces and "
775 << tgtTotalSize << " target faces"
776 << endl;
777
778 singlePatchProc_ = calcDistribution(srcPatch, tgtPatch);
779
780 if (debug)
781 {
782 Info<< "AMIInterpolation:" << nl
783 << " singlePatchProc:" << singlePatchProc_ << nl
784 << endl;
785 }
786
787 return true;
788}
789
790
792(
793 autoPtr<mapDistribute>&& srcToTgtMap,
794 autoPtr<mapDistribute>&& tgtToSrcMap,
795 labelListList&& srcAddress,
796 scalarListList&& srcWeights,
797 labelListList&& tgtAddress,
798 scalarListList&& tgtWeights
799)
800{
802
803 srcAddress_.transfer(srcAddress);
804 srcWeights_.transfer(srcWeights);
805 tgtAddress_.transfer(tgtAddress);
806 tgtWeights_.transfer(tgtWeights);
807
808 // Reset the sums of the weights
809 srcWeightsSum_.setSize(srcWeights_.size());
810 forAll(srcWeights_, facei)
811 {
812 srcWeightsSum_[facei] = sum(srcWeights_[facei]);
813 }
814
815 tgtWeightsSum_.setSize(tgtWeights_.size());
816 forAll(tgtWeights_, facei)
817 {
818 tgtWeightsSum_[facei] = sum(tgtWeights_[facei]);
819 }
820
821 srcMapPtr_ = std::move(srcToTgtMap);
822 tgtMapPtr_ = std::move(tgtToSrcMap);
823
824 upToDate_ = true;
825}
826
827
829(
830 const primitivePatch& srcPatch,
831 const primitivePatch& tgtPatch
832)
833{
834 addProfiling(ami, "AMIInterpolation::append");
835
836 // Create a new interpolation
837 auto newPtr = clone();
838 newPtr->calculate(srcPatch, tgtPatch);
839
840 // If parallel then combine the mapDistribution and re-index
841 if (distributed())
842 {
843 labelListList& srcSubMap = srcMapPtr_->subMap();
844 labelListList& srcConstructMap = srcMapPtr_->constructMap();
845
846 labelListList& tgtSubMap = tgtMapPtr_->subMap();
847 labelListList& tgtConstructMap = tgtMapPtr_->constructMap();
848
849 labelListList& newSrcSubMap = newPtr->srcMapPtr_->subMap();
850 labelListList& newSrcConstructMap = newPtr->srcMapPtr_->constructMap();
851
852 labelListList& newTgtSubMap = newPtr->tgtMapPtr_->subMap();
853 labelListList& newTgtConstructMap = newPtr->tgtMapPtr_->constructMap();
854
855 // Re-calculate the source indices
856 {
857 labelList mapMap(0), newMapMap(0);
858 forAll(srcSubMap, proci)
859 {
860 mapMap.append
861 (
863 (
864 srcConstructMap[proci].size(),
865 (mapMap.size() + newMapMap.size())
866 )
867 );
868 newMapMap.append
869 (
871 (
872 newSrcConstructMap[proci].size(),
873 (mapMap.size() + newMapMap.size())
874 )
875 );
876 }
877
878 forAll(srcSubMap, proci)
879 {
880 forAll(srcConstructMap[proci], srci)
881 {
882 srcConstructMap[proci][srci] =
883 mapMap[srcConstructMap[proci][srci]];
884 }
885 }
886
887 forAll(srcSubMap, proci)
888 {
889 forAll(newSrcConstructMap[proci], srci)
890 {
891 newSrcConstructMap[proci][srci] =
892 newMapMap[newSrcConstructMap[proci][srci]];
893 }
894 }
895
896 forAll(tgtAddress_, tgti)
897 {
898 forAll(tgtAddress_[tgti], tgtj)
899 {
900 tgtAddress_[tgti][tgtj] = mapMap[tgtAddress_[tgti][tgtj]];
901 }
902 }
903
904 forAll(newPtr->tgtAddress_, tgti)
905 {
906 forAll(newPtr->tgtAddress_[tgti], tgtj)
907 {
908 newPtr->tgtAddress_[tgti][tgtj] =
909 newMapMap[newPtr->tgtAddress_[tgti][tgtj]];
910 }
911 }
912 }
913
914 // Re-calculate the target indices
915 {
916 labelList mapMap(0), newMapMap(0);
917 forAll(srcSubMap, proci)
918 {
919 mapMap.append
920 (
922 (
923 tgtConstructMap[proci].size(),
924 (mapMap.size() + newMapMap.size())
925 )
926 );
927 newMapMap.append
928 (
930 (
931 newTgtConstructMap[proci].size(),
932 (mapMap.size() + newMapMap.size())
933 )
934 );
935 }
936
937 forAll(srcSubMap, proci)
938 {
939 forAll(tgtConstructMap[proci], tgti)
940 {
941 tgtConstructMap[proci][tgti] =
942 mapMap[tgtConstructMap[proci][tgti]];
943 }
944 }
945
946 forAll(srcSubMap, proci)
947 {
948 forAll(newTgtConstructMap[proci], tgti)
949 {
950 newTgtConstructMap[proci][tgti] =
951 newMapMap[newTgtConstructMap[proci][tgti]];
952 }
953 }
954
955 forAll(srcAddress_, srci)
956 {
957 forAll(srcAddress_[srci], srcj)
958 {
959 srcAddress_[srci][srcj] =
960 mapMap[srcAddress_[srci][srcj]];
961 }
962 }
963
964 forAll(newPtr->srcAddress_, srci)
965 {
966 forAll(newPtr->srcAddress_[srci], srcj)
967 {
968 newPtr->srcAddress_[srci][srcj] =
969 newMapMap[newPtr->srcAddress_[srci][srcj]];
970 }
971 }
972 }
973
974 // Sum the construction sizes
975 srcMapPtr_->constructSize() += newPtr->srcMapPtr_->constructSize();
976 tgtMapPtr_->constructSize() += newPtr->tgtMapPtr_->constructSize();
977
978 // Combine the maps
979 forAll(srcSubMap, proci)
980 {
981 srcSubMap[proci].append(newSrcSubMap[proci]);
982 srcConstructMap[proci].append(newSrcConstructMap[proci]);
983
984 tgtSubMap[proci].append(newTgtSubMap[proci]);
985 tgtConstructMap[proci].append(newTgtConstructMap[proci]);
986 }
987 }
988
989 // Combine new and current source data
990 forAll(srcMagSf_, srcFacei)
991 {
992 srcAddress_[srcFacei].append(newPtr->srcAddress()[srcFacei]);
993 srcWeights_[srcFacei].append(newPtr->srcWeights()[srcFacei]);
994 srcWeightsSum_[srcFacei] += newPtr->srcWeightsSum()[srcFacei];
995 }
996
997 // Combine new and current target data
998 forAll(tgtMagSf_, tgtFacei)
999 {
1000 tgtAddress_[tgtFacei].append(newPtr->tgtAddress()[tgtFacei]);
1001 tgtWeights_[tgtFacei].append(newPtr->tgtWeights()[tgtFacei]);
1002 tgtWeightsSum_[tgtFacei] += newPtr->tgtWeightsSum()[tgtFacei];
1003 }
1004}
1005
1006
1008(
1009 const bool conformal,
1010 const bool output
1011)
1012{
1013 normaliseWeights
1014 (
1015 srcMagSf_,
1016 "source",
1017 srcAddress_,
1018 srcWeights_,
1019 srcWeightsSum_,
1020 conformal,
1021 output,
1022 lowWeightCorrection_
1023 );
1024
1025 normaliseWeights
1026 (
1027 tgtMagSf_,
1028 "target",
1029 tgtAddress_,
1030 tgtWeights_,
1031 tgtWeightsSum_,
1032 conformal,
1033 output,
1034 lowWeightCorrection_
1035 );
1036}
1037
1038
1040(
1041 const primitivePatch& srcPatch,
1042 const primitivePatch& tgtPatch,
1043 const vector& n,
1044 const label tgtFacei,
1045 point& tgtPoint
1046)
1047const
1048{
1049 const pointField& srcPoints = srcPatch.points();
1050
1051 // Source face addresses that intersect target face tgtFacei
1052 const labelList& addr = tgtAddress_[tgtFacei];
1053
1054 pointHit nearest;
1055 nearest.setDistance(GREAT);
1056 label nearestFacei = -1;
1057
1058 for (const label srcFacei : addr)
1059 {
1060 const face& f = srcPatch[srcFacei];
1061
1062 pointHit ray =
1063 f.ray(tgtPoint, n, srcPoints, intersection::algorithm::VISIBLE);
1064
1065 if (ray.hit())
1066 {
1067 tgtPoint = ray.rawPoint();
1068 return srcFacei;
1069 }
1070 else if (ray.distance() < nearest.distance())
1071 {
1072
1073 nearest = ray;
1074 nearestFacei = srcFacei;
1075 }
1076 }
1077
1078 if (nearest.hit() || nearest.eligibleMiss())
1079 {
1080 tgtPoint = nearest.rawPoint();
1081 return nearestFacei;
1082 }
1083
1084 return -1;
1085}
1086
1087
1089(
1090 const primitivePatch& srcPatch,
1091 const primitivePatch& tgtPatch,
1092 const vector& n,
1093 const label srcFacei,
1094 point& srcPoint
1095)
1096const
1097{
1098 const pointField& tgtPoints = tgtPatch.points();
1099
1100 pointHit nearest;
1101 nearest.setDistance(GREAT);
1102 label nearestFacei = -1;
1103
1104 // Target face addresses that intersect source face srcFacei
1105 const labelList& addr = srcAddress_[srcFacei];
1106
1107 for (const label tgtFacei : addr)
1108 {
1109 const face& f = tgtPatch[tgtFacei];
1110
1111 pointHit ray =
1112 f.ray(srcPoint, n, tgtPoints, intersection::algorithm::VISIBLE);
1113
1114 if (ray.hit())
1115 {
1116 srcPoint = ray.rawPoint();
1117 return tgtFacei;
1118 }
1119 const pointHit near = f.nearestPoint(srcPoint, tgtPoints);
1120
1121 if (near.distance() < nearest.distance())
1122 {
1123 nearest = near;
1124 nearestFacei = tgtFacei;
1125 }
1126 }
1127 if (nearest.hit() || nearest.eligibleMiss())
1128 {
1129
1130 srcPoint = nearest.rawPoint();
1131 return nearestFacei;
1132 }
1133
1134 return -1;
1135}
1136
1137
1139{
1140 if (Pstream::parRun() && (singlePatchProc_ == -1))
1141 {
1142 Log << "Checks only valid for serial running (currently)" << endl;
1143
1144 return true;
1145 }
1146
1147 bool symmetricSrc = true;
1148
1149 Log << " Checking for missing src face in tgt lists" << nl;
1150
1151 forAll(srcAddress_, srcFacei)
1152 {
1153 const labelList& tgtIds = srcAddress_[srcFacei];
1154 for (const label tgtFacei : tgtIds)
1155 {
1156 if (!tgtAddress_[tgtFacei].found(srcFacei))
1157 {
1158 symmetricSrc = false;
1159
1160 Log << " srcFacei:" << srcFacei
1161 << " not found in tgtToSrc list for tgtFacei:"
1162 << tgtFacei << nl;
1163 }
1164 }
1165 }
1166
1167 if (symmetricSrc)
1168 {
1169 Log << " - symmetric" << endl;
1170 }
1171
1172 bool symmetricTgt = true;
1173
1174 Log << " Checking for missing tgt face in src lists" << nl;
1175
1176 forAll(tgtAddress_, tgtFacei)
1177 {
1178 const labelList& srcIds = tgtAddress_[tgtFacei];
1179 for (const label srcFacei : srcIds)
1180 {
1181 if (!srcAddress_[srcFacei].found(tgtFacei))
1182 {
1183 symmetricTgt = false;
1184
1185 Log << " tgtFacei:" << tgtFacei
1186 << " not found in srcToTgt list for srcFacei:"
1187 << srcFacei << nl;
1188 }
1189 }
1190 }
1191
1192 if (symmetricTgt)
1193 {
1194 Log << " - symmetric" << endl;
1195 }
1196
1197 return symmetricSrc && symmetricTgt;
1198}
1199
1200
1202(
1203 const primitivePatch& srcPatch,
1204 const primitivePatch& tgtPatch,
1205 const labelListList& srcAddress
1206)
1207const
1208{
1209 OFstream os("faceConnectivity" + Foam::name(Pstream::myProcNo()) + ".obj");
1210
1211 label pti = 1;
1212
1213 forAll(srcAddress, i)
1214 {
1215 const labelList& addr = srcAddress[i];
1216 const point& srcPt = srcPatch.faceCentres()[i];
1217
1218 for (const label tgtPti : addr)
1219 {
1220 const point& tgtPt = tgtPatch.faceCentres()[tgtPti];
1221
1222 meshTools::writeOBJ(os, srcPt);
1223 meshTools::writeOBJ(os, tgtPt);
1224
1225 os << "l " << pti << " " << pti + 1 << endl;
1226
1227 pti += 2;
1228 }
1229 }
1230}
1231
1232
1234{
1235 os.writeEntry("AMIMethod", type());
1236
1237 if (!requireMatch_)
1238 {
1239 os.writeEntry("requireMatch", requireMatch_);
1240 }
1241
1242 if (reverseTarget_)
1243 {
1244 os.writeEntry("reverseTarget", reverseTarget_);
1245 }
1246
1247 if (lowWeightCorrection_ > 0)
1248 {
1249 os.writeEntry("lowWeightCorrection", lowWeightCorrection_);
1250 }
1251}
1252
1253
1254// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
bool found
label n
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
labelListList srcAddress_
Addresses of target faces per source face.
scalarList tgtMagSf_
Target face areas.
const scalarListList & tgtWeights() const
Return const access to target patch weights.
void projectPointsToSurface(const searchableSurface &surf, pointField &pts) const
Project points to surface.
const labelListList & srcAddress() const
Return const access to source patch addressing.
autoPtr< mapDistribute > srcMapPtr_
Source map pointer - parallel running only.
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
bool checkSymmetricWeights(const bool log) const
autoPtr< mapDistribute > tgtMapPtr_
Target map pointer - parallel running only.
const scalarListList & srcWeights() const
Return const access to source patch weights.
static void agglomerate(const autoPtr< mapDistribute > &targetMap, const scalarList &fineSrcMagSf, const labelListList &fineSrcAddress, const scalarListList &fineSrcWeights, const labelList &sourceRestrictAddressing, const labelList &targetRestrictAddressing, scalarList &srcMagSf, labelListList &srcAddress, scalarListList &srcWeights, scalarField &srcWeightsSum, autoPtr< mapDistribute > &tgtMap)
static void normaliseWeights(const scalarList &patchAreas, const word &patchName, const labelListList &addr, scalarListList &wght, scalarField &wghtSum, const bool conformal, const bool output, const scalar lowWeightTol)
labelListList tgtAddress_
Addresses of source faces per target face.
label srcPointFace(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const vector &n, const label tgtFacei, point &tgtPoint) const
Return source patch face index of point on target patch face.
label calcDistribution(const primitivePatch &srcPatch, const primitivePatch &tgtPatch) const
Calculate if patches are on multiple processors.
static bool cacheIntersections_
void writeFaceConnectivity(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const labelListList &srcAddress) const
Write face connectivity as OBJ file.
scalarField srcWeightsSum_
Sum of weights of target faces per source face.
autoPtr< indexedOctree< treeType > > createTree(const primitivePatch &patch) const
Reset the octree for the patch face search.
scalarListList tgtWeights_
Weights of source faces per target face.
scalarListList srcWeights_
Weights of target faces per source face.
const List< scalar > & srcMagSf() const
Return const access to source patch face areas.
label tgtPointFace(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const vector &n, const label srcFacei, point &srcPoint) const
Return target patch face index of point on source patch face.
scalarList srcMagSf_
Source face areas.
scalarField tgtWeightsSum_
Sum of weights of source faces per target face.
const List< scalar > & tgtMagSf() const
Return const access to target patch face areas.
const labelListList & tgtAddress() const
Return const access to target patch addressing.
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Output to file stream, using an OSstream.
Definition: OFstream.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:54
const point_type & rawPoint() const noexcept
The point, no checks.
Definition: PointHit.H:172
bool eligibleMiss() const noexcept
Is this an eligible miss.
Definition: PointHit.H:127
scalar distance() const noexcept
Return distance to hit.
Definition: PointHit.H:139
void setDistance(const scalar d) noexcept
Set the distance.
Definition: PointHit.H:201
bool hit() const noexcept
Is there a hit.
Definition: PointHit.H:121
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
A list of faces which address into the list of points.
const Field< point_type > & points() const noexcept
Return reference to global points.
const Field< point_type > & faceCentres() const
Return face centres for patch.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void broadcast(Type &value, const label comm=UPstream::worldComm)
A List obtained as a section of another List.
Definition: SubList.H:70
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:500
label find_first() const
Locate the first bit that is set.
Definition: bitSetI.H:314
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:175
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
void reset()
Reset to defaults.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
virtual bool write()
Write the output fields.
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
const labelListList & constructMap() const noexcept
From subsetted data to new reconstructed data.
const labelListList & subMap() const noexcept
From subsetted data back to original data.
label constructSize() const noexcept
Constructed data size.
Class containing processor-to-processor mapping information.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
int myProcNo() const noexcept
Return processor number.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
Encapsulation of data needed to search on PrimitivePatches.
bool append() const noexcept
True if output format uses an append mode.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
constexpr label labelMin
Definition: label.H:60
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Type gAverage(const FieldField< Field, Type > &f)
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
labelList f(nPoints)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333