AMIInterpolation.C
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-2022 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 "AMIInterpolation.H"
27 #include "AMIMethod.H"
28 #include "meshTools.H"
29 #include "distributionMap.H"
30 #include "flipOp.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(AMIInterpolation, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
41 
43 (
44  const primitivePatch& patch,
46 )
47 {
48  tmp<scalarField> tResult(new scalarField(patch.size(), Zero));
49  scalarField& result = tResult.ref();
50 
51  const pointField& patchPoints = patch.localPoints();
52 
53  triFaceList patchFaceTris;
54 
55  forAll(result, patchFacei)
56  {
58  (
59  patch.localFaces()[patchFacei],
60  patchPoints,
61  triMode,
62  patchFaceTris
63  );
64 
65  forAll(patchFaceTris, i)
66  {
67  result[patchFacei] +=
68  patchFaceTris[i].tri(patchPoints).mag();
69  }
70  }
71 
72  return tResult;
73 }
74 
75 
76 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
77 
78 void Foam::AMIInterpolation::projectPointsToSurface
79 (
80  const searchableSurface& surf,
81  pointField& pts
82 ) const
83 {
84  if (debug)
85  {
86  Info<< "AMI: projecting points to surface" << endl;
87  }
88 
90 
91  surf.findNearest(pts, scalarField(pts.size(), great), nearInfo);
92 
93  label nMiss = 0;
94  forAll(nearInfo, i)
95  {
96  const pointIndexHit& pi = nearInfo[i];
97 
98  if (pi.hit())
99  {
100  pts[i] = pi.hitPoint();
101  }
102  else
103  {
104  pts[i] = pts[i];
105  nMiss++;
106  }
107  }
108 
109  if (nMiss > 0)
110  {
112  << "Error projecting points to surface: "
113  << nMiss << " faces could not be determined"
114  << abort(FatalError);
115  }
116 }
117 
118 
119 void Foam::AMIInterpolation::sumWeights
120 (
121  const scalarListList& wght,
122  scalarField& wghtSum
123 )
124 {
125  wghtSum.setSize(wght.size());
126  wghtSum = Zero;
127 
128  forAll(wght, facei)
129  {
130  wghtSum[facei] = sum(wght[facei]);
131  }
132 }
133 
134 
135 void Foam::AMIInterpolation::sumWeights
136 (
137  const UPtrList<scalarListList>& wghts,
138  scalarField& wghtSum
139 )
140 {
141  wghtSum.setSize(wghts[0].size());
142  wghtSum = Zero;
143 
144  forAll(wghts[0], facei)
145  {
146  forAll(wghts, wghtsi)
147  {
148  forAll(wghts[wghtsi][facei], i)
149  {
150  wghtSum[facei] += wghts[wghtsi][facei][i];
151  }
152  }
153  }
154 }
155 
156 
157 void Foam::AMIInterpolation::reportSumWeights
158 (
159  const scalarField& patchAreas,
160  const word& patchName,
161  const scalarField& wghtSum,
162  const scalar lowWeightTol
163 )
164 {
165  if (returnReduce(wghtSum.size(), sumOp<label>()) == 0)
166  {
167  return;
168  }
169 
170  label nLowWeight = 0;
171  forAll(wghtSum, facei)
172  {
173  if (wghtSum[facei] < lowWeightTol)
174  {
175  ++ nLowWeight;
176  }
177  }
178  reduce(nLowWeight, sumOp<label>());
179 
180  Info<< indent << "AMI: Patch " << patchName
181  << " sum(weights) min/max/average = " << gMin(wghtSum) << ", "
182  << gMax(wghtSum) << ", "
183  << gSum(wghtSum*patchAreas)/gSum(patchAreas) << endl;
184 
185  if (nLowWeight)
186  {
187  Info<< indent << "AMI: Patch " << patchName << " identified "
188  << nLowWeight << " faces with weights less than "
189  << lowWeightTol << endl;
190  }
191 }
192 
193 
194 void Foam::AMIInterpolation::normaliseWeights
195 (
196  scalarListList& wght,
197  const scalarField& wghtSum
198 )
199 {
200  forAll(wghtSum, facei)
201  {
202  scalarList& w = wght[facei];
203 
204  forAll(w, i)
205  {
206  w[i] /= wghtSum[facei];
207  }
208  }
209 }
210 
211 
212 void Foam::AMIInterpolation::normaliseWeights
213 (
215  const scalarField& wghtSum
216 )
217 {
218  forAll(wghtSum, facei)
219  {
220  forAll(wghts, wghtsi)
221  {
222  scalarList& w = wghts[wghtsi][facei];
223 
224  forAll(w, i)
225  {
226  w[i] /= wghtSum[facei];
227  }
228  }
229  }
230 }
231 
232 
233 void Foam::AMIInterpolation::agglomerate
234 (
235  const autoPtr<distributionMap>& targetMapPtr,
236  const scalarField& fineSrcMagSf,
237  const labelListList& fineSrcAddress,
238  const scalarListList& fineSrcWeights,
239 
240  const labelList& sourceRestrictAddressing,
241  const labelList& targetRestrictAddressing,
242 
243  scalarField& srcMagSf,
244  labelListList& srcAddress,
245  scalarListList& srcWeights,
246  scalarField& srcWeightsSum,
248 )
249 {
250  label sourceCoarseSize =
251  (
252  sourceRestrictAddressing.size()
253  ? max(sourceRestrictAddressing)+1
254  : 0
255  );
256 
257  label targetCoarseSize =
258  (
259  targetRestrictAddressing.size()
260  ? max(targetRestrictAddressing)+1
261  : 0
262  );
263 
264  // Agglomerate face areas
265  {
266  srcMagSf.setSize(sourceRestrictAddressing.size(), 0.0);
267  forAll(sourceRestrictAddressing, facei)
268  {
269  label coarseFacei = sourceRestrictAddressing[facei];
270  srcMagSf[coarseFacei] += fineSrcMagSf[facei];
271  }
272  }
273 
274 
275  // Agglomerate weights and indices
276  if (targetMapPtr.valid())
277  {
278  const distributionMap& map = targetMapPtr();
279 
280  // Get all restriction addressing.
281  labelList allRestrict(targetRestrictAddressing);
282  map.distribute(allRestrict);
283 
284  // So now we have agglomeration of the target side in
285  // allRestrict:
286  // 0..size-1 : local agglomeration (= targetRestrictAddressing)
287  // size.. : agglomeration data from other processors
288 
289  labelListList tgtSubMap(Pstream::nProcs());
290 
291  // Local subMap is just identity
292  {
293  tgtSubMap[Pstream::myProcNo()] = identity(targetCoarseSize);
294  }
295 
296  forAll(map.subMap(), proci)
297  {
298  if (proci != Pstream::myProcNo())
299  {
300  // Combine entries that point to the same coarse element. All
301  // the elements refer to local data so index into
302  // targetRestrictAddressing or allRestrict (since the same
303  // for local data).
304  const labelList& elems = map.subMap()[proci];
305  labelList& newSubMap = tgtSubMap[proci];
306  newSubMap.setSize(elems.size());
307 
308  labelList oldToNew(targetCoarseSize, -1);
309  label newI = 0;
310 
311  forAll(elems, i)
312  {
313  label fineElem = elems[i];
314  label coarseElem = allRestrict[fineElem];
315  if (oldToNew[coarseElem] == -1)
316  {
317  oldToNew[coarseElem] = newI;
318  newSubMap[newI] = coarseElem;
319  newI++;
320  }
321  }
322  newSubMap.setSize(newI);
323  }
324  }
325 
326  // Reconstruct constructMap by combining entries. Note that order
327  // of handing out indices should be the same as loop above to compact
328  // the sending map
329 
330  labelListList tgtConstructMap(Pstream::nProcs());
331  labelList tgtCompactMap;
332 
333  // Local constructMap is just identity
334  {
335  tgtConstructMap[Pstream::myProcNo()] =
336  identity(targetCoarseSize);
337 
338  tgtCompactMap = targetRestrictAddressing;
339  }
340  tgtCompactMap.setSize(map.constructSize());
341  label compactI = targetCoarseSize;
342 
343  // Compact data from other processors
344  forAll(map.constructMap(), proci)
345  {
346  if (proci != Pstream::myProcNo())
347  {
348  // Combine entries that point to the same coarse element. All
349  // elements now are remote data so we cannot use any local
350  // data here - use allRestrict instead.
351  const labelList& elems = map.constructMap()[proci];
352 
353  labelList& newConstructMap = tgtConstructMap[proci];
354  newConstructMap.setSize(elems.size());
355 
356  if (elems.size())
357  {
358  // Get the maximum target coarse size for this set of
359  // received data.
360  label remoteTargetCoarseSize = labelMin;
361  forAll(elems, i)
362  {
363  remoteTargetCoarseSize = max
364  (
365  remoteTargetCoarseSize,
366  allRestrict[elems[i]]
367  );
368  }
369  remoteTargetCoarseSize += 1;
370 
371  // Combine locally data coming from proci
372  labelList oldToNew(remoteTargetCoarseSize, -1);
373  label newI = 0;
374 
375  forAll(elems, i)
376  {
377  label fineElem = elems[i];
378  // fineElem now points to section from proci
379  label coarseElem = allRestrict[fineElem];
380  if (oldToNew[coarseElem] == -1)
381  {
382  oldToNew[coarseElem] = newI;
383  tgtCompactMap[fineElem] = compactI;
384  newConstructMap[newI] = compactI++;
385  newI++;
386  }
387  else
388  {
389  // Get compact index
390  label compactI = oldToNew[coarseElem];
391  tgtCompactMap[fineElem] = newConstructMap[compactI];
392  }
393  }
394  newConstructMap.setSize(newI);
395  }
396  }
397  }
398 
399 
400  srcAddress.setSize(sourceCoarseSize);
401  srcWeights.setSize(sourceCoarseSize);
402 
403  forAll(fineSrcAddress, facei)
404  {
405  // All the elements contributing to facei. Are slots in
406  // distributionMap'd data.
407  const labelList& elems = fineSrcAddress[facei];
408  const scalarList& weights = fineSrcWeights[facei];
409  const scalar fineArea = fineSrcMagSf[facei];
410 
411  label coarseFacei = sourceRestrictAddressing[facei];
412 
413  const scalar coarseArea = srcMagSf[coarseFacei];
414 
415  labelList& newElems = srcAddress[coarseFacei];
416  scalarList& newWeights = srcWeights[coarseFacei];
417 
418  forAll(elems, i)
419  {
420  label elemI = elems[i];
421  label coarseElemI = tgtCompactMap[elemI];
422 
423  label index = findIndex(newElems, coarseElemI);
424  if (index == -1)
425  {
426  newElems.append(coarseElemI);
427  newWeights.append(fineArea/coarseArea*weights[i]);
428  }
429  else
430  {
431  newWeights[index] += fineArea/coarseArea*weights[i];
432  }
433  }
434  }
435 
436  tgtMap.reset
437  (
438  new distributionMap
439  (
440  compactI,
441  move(tgtSubMap),
442  move(tgtConstructMap)
443  )
444  );
445  }
446  else
447  {
448  srcAddress.setSize(sourceCoarseSize);
449  srcWeights.setSize(sourceCoarseSize);
450 
451  forAll(fineSrcAddress, facei)
452  {
453  // All the elements contributing to facei. Are slots in
454  // distributionMap'd data.
455  const labelList& elems = fineSrcAddress[facei];
456  const scalarList& weights = fineSrcWeights[facei];
457  const scalar fineArea = fineSrcMagSf[facei];
458 
459  label coarseFacei = sourceRestrictAddressing[facei];
460 
461  const scalar coarseArea = srcMagSf[coarseFacei];
462 
463  labelList& newElems = srcAddress[coarseFacei];
464  scalarList& newWeights = srcWeights[coarseFacei];
465 
466  forAll(elems, i)
467  {
468  label elemI = elems[i];
469  label coarseElemI = targetRestrictAddressing[elemI];
470 
471  label index = findIndex(newElems, coarseElemI);
472  if (index == -1)
473  {
474  newElems.append(coarseElemI);
475  newWeights.append(fineArea/coarseArea*weights[i]);
476  }
477  else
478  {
479  newWeights[index] += fineArea/coarseArea*weights[i];
480  }
481  }
482  }
483  }
484 }
485 
486 
487 void Foam::AMIInterpolation::constructFromSurface
488 (
489  const primitivePatch& srcPatch,
490  const primitivePatch& tgtPatch,
491  const autoPtr<searchableSurface>& surfPtr,
492  const bool report
493 )
494 {
495  if (surfPtr.valid())
496  {
497  // Create new patches for source and target
498  pointField srcPoints = srcPatch.points();
499  primitivePatch srcPatch0
500  (
502  (
503  srcPatch,
504  srcPatch.size(),
505  0
506  ),
507  srcPoints
508  );
509 
510  if (debug)
511  {
512  OFstream os("amiSrcPoints.obj");
513  forAll(srcPoints, i)
514  {
515  meshTools::writeOBJ(os, srcPoints[i]);
516  }
517  }
518 
519  pointField tgtPoints = tgtPatch.points();
520  primitivePatch tgtPatch0
521  (
523  (
524  tgtPatch,
525  tgtPatch.size(),
526  0
527  ),
528  tgtPoints
529  );
530 
531  if (debug)
532  {
533  OFstream os("amiTgtPoints.obj");
534  forAll(tgtPoints, i)
535  {
536  meshTools::writeOBJ(os, tgtPoints[i]);
537  }
538  }
539 
540 
541  // Map source and target patches onto projection surface
542  projectPointsToSurface(surfPtr(), srcPoints);
543  projectPointsToSurface(surfPtr(), tgtPoints);
544 
545 
546  // Calculate AMI interpolation
547  update(srcPatch0, tgtPatch0, report);
548  }
549  else
550  {
551  update(srcPatch, tgtPatch, report);
552  }
553 }
554 
555 
556 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
557 
559 (
560  const primitivePatch& srcPatch,
561  const primitivePatch& tgtPatch,
563  const bool requireMatch,
564  const word& methodName,
565  const scalar lowWeightCorrection,
566  const bool reverseTarget,
567  const bool report
568 )
569 :
570  methodName_(methodName),
571  reverseTarget_(reverseTarget),
572  requireMatch_(requireMatch),
573  singlePatchProc_(-999),
574  lowWeightCorrection_(lowWeightCorrection),
575  srcAddress_(),
576  srcWeights_(),
577  srcWeightsSum_(),
578  tgtAddress_(),
579  tgtWeights_(),
580  tgtWeightsSum_(),
581  triMode_(triMode),
582  srcMapPtr_(nullptr),
583  tgtMapPtr_(nullptr)
584 {
585  update(srcPatch, tgtPatch, report);
586 }
587 
588 
590 (
591  const primitivePatch& srcPatch,
592  const primitivePatch& tgtPatch,
593  const autoPtr<searchableSurface>& surfPtr,
595  const bool requireMatch,
596  const word& methodName,
597  const scalar lowWeightCorrection,
598  const bool reverseTarget,
599  const bool report
600 )
601 :
602  methodName_(methodName),
603  reverseTarget_(reverseTarget),
604  requireMatch_(requireMatch),
605  singlePatchProc_(-999),
606  lowWeightCorrection_(lowWeightCorrection),
607  srcAddress_(),
608  srcWeights_(),
609  srcWeightsSum_(),
610  tgtAddress_(),
611  tgtWeights_(),
612  tgtWeightsSum_(),
613  triMode_(triMode),
614  srcMapPtr_(nullptr),
615  tgtMapPtr_(nullptr)
616 {
617  constructFromSurface(srcPatch, tgtPatch, surfPtr, report);
618 }
619 
620 
622 (
623  const AMIInterpolation& fineAMI,
624  const labelList& sourceRestrictAddressing,
625  const labelList& targetRestrictAddressing,
626  const bool report
627 )
628 :
629  methodName_(fineAMI.methodName_),
630  reverseTarget_(fineAMI.reverseTarget_),
631  requireMatch_(fineAMI.requireMatch_),
632  singlePatchProc_(fineAMI.singlePatchProc_),
633  lowWeightCorrection_(-1.0),
634  srcAddress_(),
635  srcWeights_(),
636  srcWeightsSum_(),
637  tgtAddress_(),
638  tgtWeights_(),
639  tgtWeightsSum_(),
640  triMode_(fineAMI.triMode_),
641  srcMapPtr_(nullptr),
642  tgtMapPtr_(nullptr)
643 {
644  label sourceCoarseSize =
645  (
646  sourceRestrictAddressing.size()
647  ? max(sourceRestrictAddressing)+1
648  : 0
649  );
650 
651  label neighbourCoarseSize =
652  (
653  targetRestrictAddressing.size()
654  ? max(targetRestrictAddressing)+1
655  : 0
656  );
657 
658  if (debug & 2)
659  {
660  Pout<< "AMI: Creating addressing and weights as agglomeration of AMI :"
661  << " source:" << fineAMI.srcAddress().size()
662  << " target:" << fineAMI.tgtAddress().size()
663  << " coarse source size:" << sourceCoarseSize
664  << " neighbour source size:" << neighbourCoarseSize
665  << endl;
666  }
667 
668  if
669  (
670  fineAMI.srcAddress().size() != sourceRestrictAddressing.size()
671  || fineAMI.tgtAddress().size() != targetRestrictAddressing.size()
672  )
673  {
675  << "Size mismatch." << nl
676  << "Source patch size:" << fineAMI.srcAddress().size() << nl
677  << "Source agglomeration size:"
678  << sourceRestrictAddressing.size() << nl
679  << "Target patch size:" << fineAMI.tgtAddress().size() << nl
680  << "Target agglomeration size:"
681  << targetRestrictAddressing.size()
682  << exit(FatalError);
683  }
684 
685  // Agglomerate addresses and weights
686  agglomerate
687  (
688  fineAMI.tgtMapPtr_,
689  fineAMI.srcMagSf(),
690  fineAMI.srcAddress(),
691  fineAMI.srcWeights(),
692 
693  sourceRestrictAddressing,
694  targetRestrictAddressing,
695 
696  srcMagSf_,
697  srcAddress_,
698  srcWeights_,
699  srcWeightsSum_,
700  tgtMapPtr_
701  );
702 
703  agglomerate
704  (
705  fineAMI.srcMapPtr_,
706  fineAMI.tgtMagSf(),
707  fineAMI.tgtAddress(),
708  fineAMI.tgtWeights(),
709 
710  targetRestrictAddressing,
711  sourceRestrictAddressing,
712 
713  tgtMagSf_,
714  tgtAddress_,
715  tgtWeights_,
716  tgtWeightsSum_,
717  srcMapPtr_
718  );
719 
720  // Weight summation and normalisation
721  sumWeights(*this);
722  if (report)
723  {
724  reportSumWeights(*this);
725  }
726  if (requireMatch_)
727  {
728  normaliseWeights(*this);
729  }
730 }
731 
732 
733 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
734 
736 {}
737 
738 
739 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
740 
742 (
743  const primitivePatch& srcPatch,
744  const primitivePatch& tgtPatch,
745  const bool report
746 )
747 {
748  label srcTotalSize = returnReduce(srcPatch.size(), sumOp<label>());
749  label tgtTotalSize = returnReduce(tgtPatch.size(), sumOp<label>());
750 
751  if (srcTotalSize == 0)
752  {
753  if (debug)
754  {
755  Info<< "AMI: no source faces present - no addressing constructed"
756  << endl;
757  }
758 
759  return;
760  }
761 
762  if (report)
763  {
764  Info<< indent
765  << "AMI: Creating addressing and weights between "
766  << srcTotalSize << " source faces and "
767  << tgtTotalSize << " target faces"
768  << endl;
769  }
770 
771  // Calculate face areas
772  srcMagSf_ = patchMagSf(srcPatch, triMode_);
773  tgtMagSf_ = patchMagSf(tgtPatch, triMode_);
774 
775  // Calculate if patches present on multiple processors
776  singlePatchProc_ = calcDistribution(srcPatch, tgtPatch);
777 
778  if (singlePatchProc_ == -1)
779  {
780  // Convert local addressing to global addressing
781  globalIndex globalSrcFaces(srcPatch.size());
782  globalIndex globalTgtFaces(tgtPatch.size());
783 
784  // Create processor map of overlapping faces. This map gets
785  // (possibly remote) faces from the tgtPatch such that they (together)
786  // cover all of the srcPatch
787  autoPtr<distributionMap> mapPtr = calcProcMap(srcPatch, tgtPatch);
788  const distributionMap& map = mapPtr();
789 
790  // Create new target patch that fully encompasses source patch
791 
792  // Faces and points
793  faceList newTgtFaces;
794  pointField newTgtPoints;
795 
796  // Original faces from tgtPatch (in globalIndexing since might be
797  // remote)
798  labelList tgtFaceIDs;
799  distributeAndMergePatches
800  (
801  map,
802  tgtPatch,
803  globalTgtFaces,
804  newTgtFaces,
805  newTgtPoints,
806  tgtFaceIDs
807  );
808 
810  newTgtPatch
811  (
813  (
814  newTgtFaces,
815  newTgtFaces.size()
816  ),
817  newTgtPoints
818  );
819  scalarField newTgtMagSf(patchMagSf(newTgtPatch, triMode_));
820 
821  // Calculate AMI interpolation
822  autoPtr<AMIMethod> AMIPtr
823  (
825  (
826  methodName_,
827  srcPatch,
828  newTgtPatch,
829  srcMagSf_,
830  newTgtMagSf,
831  triMode_,
832  reverseTarget_,
833  requireMatch_ && (lowWeightCorrection_ < 0)
834  )
835  );
836 
837  AMIPtr->calculate
838  (
839  srcAddress_,
840  srcWeights_,
841  tgtAddress_,
842  tgtWeights_
843  );
844 
845  // Now
846  // ~~~
847  // srcAddress_ : per srcPatch face a list of the newTgtPatch (not
848  // tgtPatch) faces it overlaps
849  // tgtAddress_ : per newTgtPatch (not tgtPatch) face a list of the
850  // srcPatch faces it overlaps
851 
852 
853  // Rework newTgtPatch indices into globalIndices of tgtPatch
854  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
855 
856 
857  forAll(srcAddress_, i)
858  {
859  labelList& addressing = srcAddress_[i];
860  forAll(addressing, addrI)
861  {
862  addressing[addrI] = tgtFaceIDs[addressing[addrI]];
863  }
864  }
865 
866  forAll(tgtAddress_, i)
867  {
868  labelList& addressing = tgtAddress_[i];
869  forAll(addressing, addrI)
870  {
871  addressing[addrI] = globalSrcFaces.toGlobal(addressing[addrI]);
872  }
873  }
874 
875  // Send data back to originating procs. Note that contributions
876  // from different processors get added (ListAppendEqOp)
877 
879  (
881  List<labelPair>(),
882  tgtPatch.size(),
883  map.constructMap(),
884  false, // has flip
885  map.subMap(),
886  false, // has flip
887  tgtAddress_,
889  flipOp(), // flip operation
890  labelList()
891  );
892 
894  (
896  List<labelPair>(),
897  tgtPatch.size(),
898  map.constructMap(),
899  false,
900  map.subMap(),
901  false,
902  tgtWeights_,
904  flipOp(),
905  scalarList()
906  );
907 
908  // Cache maps and reset addresses
909  List<Map<label>> cMap;
910  srcMapPtr_.reset
911  (
912  new distributionMap(globalSrcFaces, tgtAddress_, cMap)
913  );
914  tgtMapPtr_.reset
915  (
916  new distributionMap(globalTgtFaces, srcAddress_, cMap)
917  );
918 
919  if (debug)
920  {
921  writeFaceConnectivity(srcPatch, newTgtPatch, srcAddress_);
922  }
923  }
924  else
925  {
926  // Calculate AMI interpolation
927  autoPtr<AMIMethod> AMIPtr
928  (
930  (
931  methodName_,
932  srcPatch,
933  tgtPatch,
934  srcMagSf_,
935  tgtMagSf_,
936  triMode_,
937  reverseTarget_,
938  requireMatch_ && (lowWeightCorrection_ < 0)
939  )
940  );
941 
942  AMIPtr->calculate
943  (
944  srcAddress_,
945  srcWeights_,
946  tgtAddress_,
947  tgtWeights_
948  );
949  }
950 
951  // Weight summation and normalisation
952  sumWeights(*this);
953  if (report)
954  {
955  reportSumWeights(*this);
956  }
957  if (requireMatch_)
958  {
959  normaliseWeights(*this);
960  }
961 
962  if (debug)
963  {
964  Info<< "AMIIPatchToPatchnterpolation :"
965  << "Constructed addressing and weights" << nl
966  << " triMode :"
968  << " singlePatchProc:" << singlePatchProc_ << nl
969  << " srcMagSf :" << gSum(srcMagSf_) << nl
970  << " tgtMagSf :" << gSum(tgtMagSf_) << nl
971  << endl;
972  }
973 }
974 
975 
976 void Foam::AMIInterpolation::sumWeights(AMIInterpolation& AMI)
977 {
978  sumWeights(AMI.srcWeights_, AMI.srcWeightsSum_);
979  sumWeights(AMI.tgtWeights_, AMI.tgtWeightsSum_);
980 }
981 
982 
983 void Foam::AMIInterpolation::sumWeights(PtrList<AMIInterpolation>& AMIs)
984 {
985  UPtrList<scalarListList> srcWeights(AMIs.size());
986  forAll(AMIs, i)
987  {
988  srcWeights.set(i, &AMIs[i].srcWeights_);
989  }
990 
991  sumWeights(srcWeights, AMIs[0].srcWeightsSum_);
992 
993  for (label i = 1; i < AMIs.size(); ++ i)
994  {
995  AMIs[i].srcWeightsSum_ = AMIs[0].srcWeightsSum_;
996  }
997 
998  UPtrList<scalarListList> tgtWeights(AMIs.size());
999  forAll(AMIs, i)
1000  {
1001  tgtWeights.set(i, &AMIs[i].tgtWeights_);
1002  }
1003 
1004  sumWeights(tgtWeights, AMIs[0].tgtWeightsSum_);
1005 
1006  for (label i = 1; i < AMIs.size(); ++ i)
1007  {
1008  AMIs[i].tgtWeightsSum_ = AMIs[0].tgtWeightsSum_;
1009  }
1010 }
1011 
1012 
1013 void Foam::AMIInterpolation::reportSumWeights(AMIInterpolation& AMI)
1014 {
1015  reportSumWeights
1016  (
1017  AMI.srcMagSf_,
1018  "source",
1019  AMI.srcWeightsSum_,
1020  AMI.lowWeightCorrection_
1021  );
1022 
1023  reportSumWeights
1024  (
1025  AMI.tgtMagSf_,
1026  "target",
1027  AMI.tgtWeightsSum_,
1028  AMI.lowWeightCorrection_
1029  );
1030 }
1031 
1032 
1033 void Foam::AMIInterpolation::normaliseWeights(AMIInterpolation& AMI)
1034 {
1035  normaliseWeights(AMI.srcWeights_, AMI.srcWeightsSum_);
1036  normaliseWeights(AMI.tgtWeights_, AMI.tgtWeightsSum_);
1037 }
1038 
1039 
1040 void Foam::AMIInterpolation::normaliseWeights(UPtrList<AMIInterpolation>& AMIs)
1041 {
1042  UPtrList<scalarListList> srcWeights(AMIs.size());
1043  forAll(AMIs, i)
1044  {
1045  srcWeights.set(i, &AMIs[i].srcWeights_);
1046  }
1047 
1048  normaliseWeights(srcWeights, AMIs[0].srcWeightsSum_);
1049 
1050  UPtrList<scalarListList> tgtWeights(AMIs.size());
1051  forAll(AMIs, i)
1052  {
1053  tgtWeights.set(i, &AMIs[i].tgtWeights_);
1054  }
1055 
1056  normaliseWeights(tgtWeights, AMIs[0].tgtWeightsSum_);
1057 }
1058 
1059 
1062  const primitivePatch& srcPatch,
1063  const primitivePatch& tgtPatch,
1064  const vector& n,
1065  const label tgtFacei,
1066  point& tgtPoint
1067 )
1068 const
1069 {
1070  const pointField& srcPoints = srcPatch.points();
1071 
1072  // Source face addresses that intersect target face tgtFacei
1073  const labelList& addr = tgtAddress_[tgtFacei];
1074 
1075  pointHit nearest;
1076  label nearestFacei = -1;
1077 
1078  forAll(addr, i)
1079  {
1080  const label srcFacei = addr[i];
1081  const face& f = srcPatch[srcFacei];
1082 
1083  const pointHit ray =
1084  f.ray(tgtPoint, n, srcPoints, intersection::algorithm::visible);
1085 
1086  if (ray.hit())
1087  {
1088  tgtPoint = ray.rawPoint();
1089  return srcFacei;
1090  }
1091 
1092  const pointHit near = f.nearestPoint(tgtPoint, srcPoints);
1093 
1094  if (near.distance() < nearest.distance())
1095  {
1096  nearest = near;
1097  nearestFacei = srcFacei;
1098  }
1099  }
1100 
1101  if (nearest.hit() || nearest.eligibleMiss())
1102  {
1103  tgtPoint = nearest.rawPoint();
1104  return nearestFacei;
1105  }
1106 
1107  return -1;
1108 }
1109 
1110 
1113  const primitivePatch& srcPatch,
1114  const primitivePatch& tgtPatch,
1115  const vector& n,
1116  const label srcFacei,
1117  point& srcPoint
1118 )
1119 const
1120 {
1121  const pointField& tgtPoints = tgtPatch.points();
1122 
1123  pointHit nearest;
1124  label nearestFacei = -1;
1125 
1126  // Target face addresses that intersect source face srcFacei
1127  const labelList& addr = srcAddress_[srcFacei];
1128 
1129  forAll(addr, i)
1130  {
1131  const label tgtFacei = addr[i];
1132  const face& f = tgtPatch[tgtFacei];
1133 
1134  const pointHit ray =
1135  f.ray(srcPoint, n, tgtPoints, intersection::algorithm::visible);
1136 
1137  if (ray.hit())
1138  {
1139  srcPoint = ray.rawPoint();
1140  return tgtFacei;
1141  }
1142 
1143  const pointHit near = f.nearestPoint(srcPoint, tgtPoints);
1144 
1145  if (near.distance() < nearest.distance())
1146  {
1147  nearest = near;
1148  nearestFacei = tgtFacei;
1149  }
1150  }
1151 
1152  if (nearest.hit() || nearest.eligibleMiss())
1153  {
1154  srcPoint = nearest.rawPoint();
1155  return nearestFacei;
1156  }
1157 
1158  return -1;
1159 }
1160 
1161 
1164  const primitivePatch& srcPatch,
1165  const primitivePatch& tgtPatch,
1166  const labelListList& srcAddress
1167 )
1168 const
1169 {
1170  OFstream os("faceConnectivity" + Foam::name(Pstream::myProcNo()) + ".obj");
1171 
1172  label ptI = 1;
1173 
1174  forAll(srcAddress, i)
1175  {
1176  const labelList& addr = srcAddress[i];
1177  const point& srcPt = srcPatch.faceCentres()[i];
1178 
1179  forAll(addr, j)
1180  {
1181  label tgtPtI = addr[j];
1182  const point& tgtPt = tgtPatch.faceCentres()[tgtPtI];
1183 
1184  meshTools::writeOBJ(os, srcPt);
1185  meshTools::writeOBJ(os, tgtPt);
1186 
1187  os << "l " << ptI << " " << ptI + 1 << endl;
1188 
1189  ptI += 2;
1190  }
1191  }
1192 }
1193 
1194 
1195 // ************************************************************************* //
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &, const negateOp &negOp, const int tag=UPstream::msgType())
Distribute data. Note:schedule only used for.
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
void update(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const bool report)
Update addressing and weights.
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
const scalarListList & tgtWeights() const
Return const access to target patch weights.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const Field< PointType > & faceCentres() const
Return face centres for patch.
Type gMin(const FieldField< Field, Type > &f)
void writeFaceConnectivity(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const labelListList &srcAddress) const
Write face connectivity as OBJ file.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Output to file stream.
Definition: OFstream.H:82
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.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const labelListList & tgtAddress() const
Return const access to target patch addressing.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
pointHit ray(const point &p, const vector &n, const pointField &, const intersection::algorithm alg=intersection::algorithm::fullRay, const intersection::direction dir=intersection::direction::vector) const
Return potential intersection with face with a ray starting.
static tmp< scalarField > patchMagSf(const primitivePatch &patch, const faceAreaIntersect::triangulationMode triMode)
Calculate the patch face magnitudes for the given tri-mode.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
virtual void calculate(labelListList &srcAddress, scalarListList &srcWeights, labelListList &tgtAddress, scalarListList &tgtWeights, label srcFacei=-1, label tgtFacei=-1)=0
Update addressing and weights.
const labelListList & subMap() const
From subsetted data back to original data.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
static const NamedEnum< triangulationMode, 2 > triangulationModeNames_
static void triangulate(const face &f, const pointField &points, const triangulationMode &triMode, triFaceList &faceTris)
Triangulate a face using the given triangulation mode.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static autoPtr< AMIMethod > New(const word &methodName, const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const scalarField &srcMagSf, const scalarField &tgtMagSf, const faceAreaIntersect::triangulationMode &triMode, const bool reverseTarget, const bool requireMatch)
Selector.
Definition: AMIMethodNew.C:31
const scalarField & srcMagSf() const
Return const access to source patch face areas.
A list of faces which address into the list of points.
const Point & hitPoint() const
Return hit point.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A List obtained as a section of another List.
Definition: SubList.H:53
Type gSum(const FieldField< Field, Type > &f)
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
bool hit() const
Is there a hit.
const Field< PointType > & points() const
Return reference to global points.
bool hit() const
Is there a hit.
Definition: PointHit.H:120
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
virtual ~AMIInterpolation()
Destructor.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
errorManip< error > abort(error &err)
Definition: errorManip.H:131
label constructSize() const
Constructed data size.
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
const scalarListList & srcWeights() const
Return const access to source patch weights.
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
Class containing processor-to-processor mapping information.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
AMIInterpolation(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const faceAreaIntersect::triangulationMode &triMode, const bool requireMatch=true, const word &methodName=faceAreaWeightAMI::typeName, const scalar lowWeightCorrection=-1, const bool reverseTarget=false, const bool report=true)
Construct from components.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
A class for managing temporary objects.
Definition: PtrList.H:53
bool eligibleMiss() const
Is this an eligible miss.
Definition: PointHit.H:164
const labelListList & srcAddress() const
Return const access to source patch addressing.
static const label labelMin
Definition: label.H:61
const scalarField & tgtMagSf() const
Return const access to target patch face areas.
List operator to append one list onto another.
Definition: ListOps.H:306
Namespace for OpenFOAM.