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-2019 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 "mapDistribute.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 interpolationMethod& im
45 )
46 {
47  word method = "unknown-interpolationMethod";
48 
49  switch (im)
50  {
51  case imDirect:
52  {
53  method = "directAMI";
54  break;
55  }
56  case imMapNearest:
57  {
58  method = "mapNearestAMI";
59  break;
60  }
61  case imFaceAreaWeight:
62  {
63  method = "faceAreaWeightAMI";
64  break;
65  }
66  case imPartialFaceAreaWeight:
67  {
68  method = "partialFaceAreaWeightAMI";
69  break;
70  }
71  case imSweptFaceAreaWeight:
72  {
73  method = "sweptFaceAreaWeightAMI";
74  break;
75  }
76  default:
77  {
79  << "Unhandled interpolationMethod enumeration " << method
80  << abort(FatalError);
81  }
82  }
83 
84  return method;
85 }
86 
87 
90 {
91  interpolationMethod method = imDirect;
92 
93  wordList methods
94  (
96  (
97  "("
98  "directAMI "
99  "mapNearestAMI "
100  "faceAreaWeightAMI "
101  "partialFaceAreaWeightAMI "
102  "sweptFaceAreaWeightAMI"
103  ")"
104  )()
105  );
106 
107  if (im == "directAMI")
108  {
109  method = imDirect;
110  }
111  else if (im == "mapNearestAMI")
112  {
113  method = imMapNearest;
114  }
115  else if (im == "faceAreaWeightAMI")
116  {
117  method = imFaceAreaWeight;
118  }
119  else if (im == "partialFaceAreaWeightAMI")
120  {
121  method = imPartialFaceAreaWeight;
122  }
123  else if (im == "sweptFaceAreaWeightAMI")
124  {
125  method = imSweptFaceAreaWeight;
126  }
127  else
128  {
130  << "Invalid interpolationMethod " << im
131  << ". Valid methods are:" << methods
132  << exit(FatalError);
133  }
134 
135  return method;
136 }
137 
138 
140 (
141  const primitivePatch& patch,
143 )
144 {
145  tmp<scalarField> tResult(new scalarField(patch.size(), Zero));
146  scalarField& result = tResult.ref();
147 
148  const pointField& patchPoints = patch.localPoints();
149 
150  faceList patchFaceTris;
151 
152  forAll(result, patchFacei)
153  {
155  (
156  patch.localFaces()[patchFacei],
157  patchPoints,
158  triMode,
159  patchFaceTris
160  );
161 
162  forAll(patchFaceTris, i)
163  {
164  result[patchFacei] +=
166  (
167  patchPoints[patchFaceTris[i][0]],
168  patchPoints[patchFaceTris[i][1]],
169  patchPoints[patchFaceTris[i][2]]
170  ).mag();
171  }
172  }
173 
174  return tResult;
175 }
176 
177 
178 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
179 
180 void Foam::AMIInterpolation::projectPointsToSurface
181 (
182  const searchableSurface& surf,
183  pointField& pts
184 ) const
185 {
186  if (debug)
187  {
188  Info<< "AMI: projecting points to surface" << endl;
189  }
190 
191  List<pointIndexHit> nearInfo;
192 
193  surf.findNearest(pts, scalarField(pts.size(), great), nearInfo);
194 
195  label nMiss = 0;
196  forAll(nearInfo, i)
197  {
198  const pointIndexHit& pi = nearInfo[i];
199 
200  if (pi.hit())
201  {
202  pts[i] = pi.hitPoint();
203  }
204  else
205  {
206  pts[i] = pts[i];
207  nMiss++;
208  }
209  }
210 
211  if (nMiss > 0)
212  {
214  << "Error projecting points to surface: "
215  << nMiss << " faces could not be determined"
216  << abort(FatalError);
217  }
218 }
219 
220 
221 void Foam::AMIInterpolation::sumWeights
222 (
223  const scalarListList& wght,
224  scalarField& wghtSum
225 )
226 {
227  wghtSum.setSize(wght.size());
228  wghtSum = Zero;
229 
230  forAll(wght, facei)
231  {
232  wghtSum[facei] = sum(wght[facei]);
233  }
234 }
235 
236 
237 void Foam::AMIInterpolation::sumWeights
238 (
239  const UPtrList<scalarListList>& wghts,
240  scalarField& wghtSum
241 )
242 {
243  wghtSum.setSize(wghts[0].size());
244  wghtSum = Zero;
245 
246  forAll(wghts[0], facei)
247  {
248  forAll(wghts, wghtsi)
249  {
250  forAll(wghts[wghtsi][facei], i)
251  {
252  wghtSum[facei] += wghts[wghtsi][facei][i];
253  }
254  }
255  }
256 }
257 
258 
259 void Foam::AMIInterpolation::reportSumWeights
260 (
261  const scalarField& patchAreas,
262  const word& patchName,
263  const scalarField& wghtSum,
264  const scalar lowWeightTol
265 )
266 {
267  if (returnReduce(wghtSum.size(), sumOp<label>()) == 0)
268  {
269  return;
270  }
271 
272  label nLowWeight = 0;
273  forAll(wghtSum, facei)
274  {
275  if (wghtSum[facei] < lowWeightTol)
276  {
277  ++ nLowWeight;
278  }
279  }
280  reduce(nLowWeight, sumOp<label>());
281 
282  Info<< indent << "AMI: Patch " << patchName
283  << " sum(weights) min/max/average = " << gMin(wghtSum) << ", "
284  << gMax(wghtSum) << ", "
285  << gSum(wghtSum*patchAreas)/gSum(patchAreas) << endl;
286 
287  if (nLowWeight)
288  {
289  Info<< indent << "AMI: Patch " << patchName << " identified "
290  << nLowWeight << " faces with weights less than "
291  << lowWeightTol << endl;
292  }
293 }
294 
295 
296 void Foam::AMIInterpolation::normaliseWeights
297 (
298  scalarListList& wght,
299  const scalarField& wghtSum
300 )
301 {
302  forAll(wghtSum, facei)
303  {
304  scalarList& w = wght[facei];
305 
306  forAll(w, i)
307  {
308  w[i] /= wghtSum[facei];
309  }
310  }
311 }
312 
313 
314 void Foam::AMIInterpolation::normaliseWeights
315 (
317  const scalarField& wghtSum
318 )
319 {
320  forAll(wghtSum, facei)
321  {
322  forAll(wghts, wghtsi)
323  {
324  scalarList& w = wghts[wghtsi][facei];
325 
326  forAll(w, i)
327  {
328  w[i] /= wghtSum[facei];
329  }
330  }
331  }
332 }
333 
334 
335 void Foam::AMIInterpolation::agglomerate
336 (
337  const autoPtr<mapDistribute>& targetMapPtr,
338  const scalarField& fineSrcMagSf,
339  const labelListList& fineSrcAddress,
340  const scalarListList& fineSrcWeights,
341 
342  const labelList& sourceRestrictAddressing,
343  const labelList& targetRestrictAddressing,
344 
345  scalarField& srcMagSf,
346  labelListList& srcAddress,
347  scalarListList& srcWeights,
348  scalarField& srcWeightsSum,
349  autoPtr<mapDistribute>& tgtMap
350 )
351 {
352  label sourceCoarseSize =
353  (
354  sourceRestrictAddressing.size()
355  ? max(sourceRestrictAddressing)+1
356  : 0
357  );
358 
359  label targetCoarseSize =
360  (
361  targetRestrictAddressing.size()
362  ? max(targetRestrictAddressing)+1
363  : 0
364  );
365 
366  // Agglomerate face areas
367  {
368  srcMagSf.setSize(sourceRestrictAddressing.size(), 0.0);
369  forAll(sourceRestrictAddressing, facei)
370  {
371  label coarseFacei = sourceRestrictAddressing[facei];
372  srcMagSf[coarseFacei] += fineSrcMagSf[facei];
373  }
374  }
375 
376 
377  // Agglomerate weights and indices
378  if (targetMapPtr.valid())
379  {
380  const mapDistribute& map = targetMapPtr();
381 
382  // Get all restriction addressing.
383  labelList allRestrict(targetRestrictAddressing);
384  map.distribute(allRestrict);
385 
386  // So now we have agglomeration of the target side in
387  // allRestrict:
388  // 0..size-1 : local agglomeration (= targetRestrictAddressing)
389  // size.. : agglomeration data from other processors
390 
391  labelListList tgtSubMap(Pstream::nProcs());
392 
393  // Local subMap is just identity
394  {
395  tgtSubMap[Pstream::myProcNo()] = identity(targetCoarseSize);
396  }
397 
398  forAll(map.subMap(), proci)
399  {
400  if (proci != Pstream::myProcNo())
401  {
402  // Combine entries that point to the same coarse element. All
403  // the elements refer to local data so index into
404  // targetRestrictAddressing or allRestrict (since the same
405  // for local data).
406  const labelList& elems = map.subMap()[proci];
407  labelList& newSubMap = tgtSubMap[proci];
408  newSubMap.setSize(elems.size());
409 
410  labelList oldToNew(targetCoarseSize, -1);
411  label newI = 0;
412 
413  forAll(elems, i)
414  {
415  label fineElem = elems[i];
416  label coarseElem = allRestrict[fineElem];
417  if (oldToNew[coarseElem] == -1)
418  {
419  oldToNew[coarseElem] = newI;
420  newSubMap[newI] = coarseElem;
421  newI++;
422  }
423  }
424  newSubMap.setSize(newI);
425  }
426  }
427 
428  // Reconstruct constructMap by combining entries. Note that order
429  // of handing out indices should be the same as loop above to compact
430  // the sending map
431 
432  labelListList tgtConstructMap(Pstream::nProcs());
433  labelList tgtCompactMap;
434 
435  // Local constructMap is just identity
436  {
437  tgtConstructMap[Pstream::myProcNo()] =
438  identity(targetCoarseSize);
439 
440  tgtCompactMap = targetRestrictAddressing;
441  }
442  tgtCompactMap.setSize(map.constructSize());
443  label compactI = targetCoarseSize;
444 
445  // Compact data from other processors
446  forAll(map.constructMap(), proci)
447  {
448  if (proci != Pstream::myProcNo())
449  {
450  // Combine entries that point to the same coarse element. All
451  // elements now are remote data so we cannot use any local
452  // data here - use allRestrict instead.
453  const labelList& elems = map.constructMap()[proci];
454 
455  labelList& newConstructMap = tgtConstructMap[proci];
456  newConstructMap.setSize(elems.size());
457 
458  if (elems.size())
459  {
460  // Get the maximum target coarse size for this set of
461  // received data.
462  label remoteTargetCoarseSize = labelMin;
463  forAll(elems, i)
464  {
465  remoteTargetCoarseSize = max
466  (
467  remoteTargetCoarseSize,
468  allRestrict[elems[i]]
469  );
470  }
471  remoteTargetCoarseSize += 1;
472 
473  // Combine locally data coming from proci
474  labelList oldToNew(remoteTargetCoarseSize, -1);
475  label newI = 0;
476 
477  forAll(elems, i)
478  {
479  label fineElem = elems[i];
480  // fineElem now points to section from proci
481  label coarseElem = allRestrict[fineElem];
482  if (oldToNew[coarseElem] == -1)
483  {
484  oldToNew[coarseElem] = newI;
485  tgtCompactMap[fineElem] = compactI;
486  newConstructMap[newI] = compactI++;
487  newI++;
488  }
489  else
490  {
491  // Get compact index
492  label compactI = oldToNew[coarseElem];
493  tgtCompactMap[fineElem] = newConstructMap[compactI];
494  }
495  }
496  newConstructMap.setSize(newI);
497  }
498  }
499  }
500 
501 
502  srcAddress.setSize(sourceCoarseSize);
503  srcWeights.setSize(sourceCoarseSize);
504 
505  forAll(fineSrcAddress, facei)
506  {
507  // All the elements contributing to facei. Are slots in
508  // mapDistribute'd data.
509  const labelList& elems = fineSrcAddress[facei];
510  const scalarList& weights = fineSrcWeights[facei];
511  const scalar fineArea = fineSrcMagSf[facei];
512 
513  label coarseFacei = sourceRestrictAddressing[facei];
514 
515  const scalar coarseArea = srcMagSf[coarseFacei];
516 
517  labelList& newElems = srcAddress[coarseFacei];
518  scalarList& newWeights = srcWeights[coarseFacei];
519 
520  forAll(elems, i)
521  {
522  label elemI = elems[i];
523  label coarseElemI = tgtCompactMap[elemI];
524 
525  label index = findIndex(newElems, coarseElemI);
526  if (index == -1)
527  {
528  newElems.append(coarseElemI);
529  newWeights.append(fineArea/coarseArea*weights[i]);
530  }
531  else
532  {
533  newWeights[index] += fineArea/coarseArea*weights[i];
534  }
535  }
536  }
537 
538  tgtMap.reset
539  (
540  new mapDistribute
541  (
542  compactI,
543  move(tgtSubMap),
544  move(tgtConstructMap)
545  )
546  );
547  }
548  else
549  {
550  srcAddress.setSize(sourceCoarseSize);
551  srcWeights.setSize(sourceCoarseSize);
552 
553  forAll(fineSrcAddress, facei)
554  {
555  // All the elements contributing to facei. Are slots in
556  // mapDistribute'd data.
557  const labelList& elems = fineSrcAddress[facei];
558  const scalarList& weights = fineSrcWeights[facei];
559  const scalar fineArea = fineSrcMagSf[facei];
560 
561  label coarseFacei = sourceRestrictAddressing[facei];
562 
563  const scalar coarseArea = srcMagSf[coarseFacei];
564 
565  labelList& newElems = srcAddress[coarseFacei];
566  scalarList& newWeights = srcWeights[coarseFacei];
567 
568  forAll(elems, i)
569  {
570  label elemI = elems[i];
571  label coarseElemI = targetRestrictAddressing[elemI];
572 
573  label index = findIndex(newElems, coarseElemI);
574  if (index == -1)
575  {
576  newElems.append(coarseElemI);
577  newWeights.append(fineArea/coarseArea*weights[i]);
578  }
579  else
580  {
581  newWeights[index] += fineArea/coarseArea*weights[i];
582  }
583  }
584  }
585  }
586 }
587 
588 
589 void Foam::AMIInterpolation::constructFromSurface
590 (
591  const primitivePatch& srcPatch,
592  const primitivePatch& tgtPatch,
593  const autoPtr<searchableSurface>& surfPtr,
594  const bool report
595 )
596 {
597  if (surfPtr.valid())
598  {
599  // Create new patches for source and target
600  pointField srcPoints = srcPatch.points();
601  primitivePatch srcPatch0
602  (
604  (
605  srcPatch,
606  srcPatch.size(),
607  0
608  ),
609  srcPoints
610  );
611 
612  if (debug)
613  {
614  OFstream os("amiSrcPoints.obj");
615  forAll(srcPoints, i)
616  {
617  meshTools::writeOBJ(os, srcPoints[i]);
618  }
619  }
620 
621  pointField tgtPoints = tgtPatch.points();
622  primitivePatch tgtPatch0
623  (
625  (
626  tgtPatch,
627  tgtPatch.size(),
628  0
629  ),
630  tgtPoints
631  );
632 
633  if (debug)
634  {
635  OFstream os("amiTgtPoints.obj");
636  forAll(tgtPoints, i)
637  {
638  meshTools::writeOBJ(os, tgtPoints[i]);
639  }
640  }
641 
642 
643  // Map source and target patches onto projection surface
644  projectPointsToSurface(surfPtr(), srcPoints);
645  projectPointsToSurface(surfPtr(), tgtPoints);
646 
647 
648  // Calculate AMI interpolation
649  update(srcPatch0, tgtPatch0, report);
650  }
651  else
652  {
653  update(srcPatch, tgtPatch, report);
654  }
655 }
656 
657 
658 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
659 
661 (
662  const primitivePatch& srcPatch,
663  const primitivePatch& tgtPatch,
665  const bool requireMatch,
666  const interpolationMethod& method,
667  const scalar lowWeightCorrection,
668  const bool reverseTarget,
669  const bool report
670 )
671 :
672  methodName_(interpolationMethodToWord(method)),
673  reverseTarget_(reverseTarget),
674  requireMatch_(requireMatch),
675  singlePatchProc_(-999),
676  lowWeightCorrection_(lowWeightCorrection),
677  srcAddress_(),
678  srcWeights_(),
679  srcWeightsSum_(),
680  tgtAddress_(),
681  tgtWeights_(),
682  tgtWeightsSum_(),
683  triMode_(triMode),
684  srcMapPtr_(nullptr),
685  tgtMapPtr_(nullptr)
686 {
687  update(srcPatch, tgtPatch, report);
688 }
689 
690 
692 (
693  const primitivePatch& srcPatch,
694  const primitivePatch& tgtPatch,
696  const bool requireMatch,
697  const word& methodName,
698  const scalar lowWeightCorrection,
699  const bool reverseTarget,
700  const bool report
701 )
702 :
703  methodName_(methodName),
704  reverseTarget_(reverseTarget),
705  requireMatch_(requireMatch),
706  singlePatchProc_(-999),
707  lowWeightCorrection_(lowWeightCorrection),
708  srcAddress_(),
709  srcWeights_(),
710  srcWeightsSum_(),
711  tgtAddress_(),
712  tgtWeights_(),
713  tgtWeightsSum_(),
714  triMode_(triMode),
715  srcMapPtr_(nullptr),
716  tgtMapPtr_(nullptr)
717 {
718  update(srcPatch, tgtPatch, report);
719 }
720 
721 
723 (
724  const primitivePatch& srcPatch,
725  const primitivePatch& tgtPatch,
726  const autoPtr<searchableSurface>& surfPtr,
728  const bool requireMatch,
729  const interpolationMethod& method,
730  const scalar lowWeightCorrection,
731  const bool reverseTarget,
732  const bool report
733 )
734 :
735  methodName_(interpolationMethodToWord(method)),
736  reverseTarget_(reverseTarget),
737  requireMatch_(requireMatch),
738  singlePatchProc_(-999),
739  lowWeightCorrection_(lowWeightCorrection),
740  srcAddress_(),
741  srcWeights_(),
742  srcWeightsSum_(),
743  tgtAddress_(),
744  tgtWeights_(),
745  tgtWeightsSum_(),
746  triMode_(triMode),
747  srcMapPtr_(nullptr),
748  tgtMapPtr_(nullptr)
749 {
750  constructFromSurface(srcPatch, tgtPatch, surfPtr, report);
751 }
752 
753 
755 (
756  const primitivePatch& srcPatch,
757  const primitivePatch& tgtPatch,
758  const autoPtr<searchableSurface>& surfPtr,
760  const bool requireMatch,
761  const word& methodName,
762  const scalar lowWeightCorrection,
763  const bool reverseTarget,
764  const bool report
765 )
766 :
767  methodName_(methodName),
768  reverseTarget_(reverseTarget),
769  requireMatch_(requireMatch),
770  singlePatchProc_(-999),
771  lowWeightCorrection_(lowWeightCorrection),
772  srcAddress_(),
773  srcWeights_(),
774  srcWeightsSum_(),
775  tgtAddress_(),
776  tgtWeights_(),
777  tgtWeightsSum_(),
778  triMode_(triMode),
779  srcMapPtr_(nullptr),
780  tgtMapPtr_(nullptr)
781 {
782  constructFromSurface(srcPatch, tgtPatch, surfPtr, report);
783 }
784 
785 
787 (
788  const AMIInterpolation& fineAMI,
789  const labelList& sourceRestrictAddressing,
790  const labelList& targetRestrictAddressing,
791  const bool report
792 )
793 :
794  methodName_(fineAMI.methodName_),
795  reverseTarget_(fineAMI.reverseTarget_),
796  requireMatch_(fineAMI.requireMatch_),
797  singlePatchProc_(fineAMI.singlePatchProc_),
798  lowWeightCorrection_(-1.0),
799  srcAddress_(),
800  srcWeights_(),
801  srcWeightsSum_(),
802  tgtAddress_(),
803  tgtWeights_(),
804  tgtWeightsSum_(),
805  triMode_(fineAMI.triMode_),
806  srcMapPtr_(nullptr),
807  tgtMapPtr_(nullptr)
808 {
809  label sourceCoarseSize =
810  (
811  sourceRestrictAddressing.size()
812  ? max(sourceRestrictAddressing)+1
813  : 0
814  );
815 
816  label neighbourCoarseSize =
817  (
818  targetRestrictAddressing.size()
819  ? max(targetRestrictAddressing)+1
820  : 0
821  );
822 
823  if (debug & 2)
824  {
825  Pout<< "AMI: Creating addressing and weights as agglomeration of AMI :"
826  << " source:" << fineAMI.srcAddress().size()
827  << " target:" << fineAMI.tgtAddress().size()
828  << " coarse source size:" << sourceCoarseSize
829  << " neighbour source size:" << neighbourCoarseSize
830  << endl;
831  }
832 
833  if
834  (
835  fineAMI.srcAddress().size() != sourceRestrictAddressing.size()
836  || fineAMI.tgtAddress().size() != targetRestrictAddressing.size()
837  )
838  {
840  << "Size mismatch." << nl
841  << "Source patch size:" << fineAMI.srcAddress().size() << nl
842  << "Source agglomeration size:"
843  << sourceRestrictAddressing.size() << nl
844  << "Target patch size:" << fineAMI.tgtAddress().size() << nl
845  << "Target agglomeration size:"
846  << targetRestrictAddressing.size()
847  << exit(FatalError);
848  }
849 
850  // Agglomerate addresses and weights
851  agglomerate
852  (
853  fineAMI.tgtMapPtr_,
854  fineAMI.srcMagSf(),
855  fineAMI.srcAddress(),
856  fineAMI.srcWeights(),
857 
858  sourceRestrictAddressing,
859  targetRestrictAddressing,
860 
861  srcMagSf_,
862  srcAddress_,
863  srcWeights_,
864  srcWeightsSum_,
865  tgtMapPtr_
866  );
867 
868  agglomerate
869  (
870  fineAMI.srcMapPtr_,
871  fineAMI.tgtMagSf(),
872  fineAMI.tgtAddress(),
873  fineAMI.tgtWeights(),
874 
875  targetRestrictAddressing,
876  sourceRestrictAddressing,
877 
878  tgtMagSf_,
879  tgtAddress_,
880  tgtWeights_,
881  tgtWeightsSum_,
882  srcMapPtr_
883  );
884 
885  // Weight summation and normalisation
886  sumWeights(*this);
887  if (report)
888  {
889  reportSumWeights(*this);
890  }
891  if (requireMatch_)
892  {
893  normaliseWeights(*this);
894  }
895 }
896 
897 
898 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
899 
901 {}
902 
903 
904 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
905 
907 (
908  const primitivePatch& srcPatch,
909  const primitivePatch& tgtPatch,
910  const bool report
911 )
912 {
913  label srcTotalSize = returnReduce(srcPatch.size(), sumOp<label>());
914  label tgtTotalSize = returnReduce(tgtPatch.size(), sumOp<label>());
915 
916  if (srcTotalSize == 0)
917  {
918  if (debug)
919  {
920  Info<< "AMI: no source faces present - no addressing constructed"
921  << endl;
922  }
923 
924  return;
925  }
926 
927  if (report)
928  {
929  Info<< indent
930  << "AMI: Creating addressing and weights between "
931  << srcTotalSize << " source faces and "
932  << tgtTotalSize << " target faces"
933  << endl;
934  }
935 
936  // Calculate face areas
937  srcMagSf_ = patchMagSf(srcPatch, triMode_);
938  tgtMagSf_ = patchMagSf(tgtPatch, triMode_);
939 
940  // Calculate if patches present on multiple processors
941  singlePatchProc_ = calcDistribution(srcPatch, tgtPatch);
942 
943  if (singlePatchProc_ == -1)
944  {
945  // Convert local addressing to global addressing
946  globalIndex globalSrcFaces(srcPatch.size());
947  globalIndex globalTgtFaces(tgtPatch.size());
948 
949  // Create processor map of overlapping faces. This map gets
950  // (possibly remote) faces from the tgtPatch such that they (together)
951  // cover all of the srcPatch
952  autoPtr<mapDistribute> mapPtr = calcProcMap(srcPatch, tgtPatch);
953  const mapDistribute& map = mapPtr();
954 
955  // Create new target patch that fully encompasses source patch
956 
957  // Faces and points
958  faceList newTgtFaces;
959  pointField newTgtPoints;
960 
961  // Original faces from tgtPatch (in globalIndexing since might be
962  // remote)
963  labelList tgtFaceIDs;
964  distributeAndMergePatches
965  (
966  map,
967  tgtPatch,
968  globalTgtFaces,
969  newTgtFaces,
970  newTgtPoints,
971  tgtFaceIDs
972  );
973 
975  newTgtPatch
976  (
978  (
979  newTgtFaces,
980  newTgtFaces.size()
981  ),
982  newTgtPoints
983  );
984  scalarField newTgtMagSf(patchMagSf(newTgtPatch, triMode_));
985 
986  // Calculate AMI interpolation
987  autoPtr<AMIMethod> AMIPtr
988  (
990  (
991  methodName_,
992  srcPatch,
993  newTgtPatch,
994  srcMagSf_,
995  newTgtMagSf,
996  triMode_,
997  reverseTarget_,
998  requireMatch_ && (lowWeightCorrection_ < 0)
999  )
1000  );
1001 
1002  AMIPtr->calculate
1003  (
1004  srcAddress_,
1005  srcWeights_,
1006  tgtAddress_,
1007  tgtWeights_
1008  );
1009 
1010  // Now
1011  // ~~~
1012  // srcAddress_ : per srcPatch face a list of the newTgtPatch (not
1013  // tgtPatch) faces it overlaps
1014  // tgtAddress_ : per newTgtPatch (not tgtPatch) face a list of the
1015  // srcPatch faces it overlaps
1016 
1017 
1018  // Rework newTgtPatch indices into globalIndices of tgtPatch
1019  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1020 
1021 
1022  forAll(srcAddress_, i)
1023  {
1024  labelList& addressing = srcAddress_[i];
1025  forAll(addressing, addrI)
1026  {
1027  addressing[addrI] = tgtFaceIDs[addressing[addrI]];
1028  }
1029  }
1030 
1031  forAll(tgtAddress_, i)
1032  {
1033  labelList& addressing = tgtAddress_[i];
1034  forAll(addressing, addrI)
1035  {
1036  addressing[addrI] = globalSrcFaces.toGlobal(addressing[addrI]);
1037  }
1038  }
1039 
1040  // Send data back to originating procs. Note that contributions
1041  // from different processors get added (ListAppendEqOp)
1042 
1044  (
1046  List<labelPair>(),
1047  tgtPatch.size(),
1048  map.constructMap(),
1049  false, // has flip
1050  map.subMap(),
1051  false, // has flip
1052  tgtAddress_,
1054  flipOp(), // flip operation
1055  labelList()
1056  );
1057 
1059  (
1061  List<labelPair>(),
1062  tgtPatch.size(),
1063  map.constructMap(),
1064  false,
1065  map.subMap(),
1066  false,
1067  tgtWeights_,
1069  flipOp(),
1070  scalarList()
1071  );
1072 
1073  // Cache maps and reset addresses
1074  List<Map<label>> cMap;
1075  srcMapPtr_.reset(new mapDistribute(globalSrcFaces, tgtAddress_, cMap));
1076  tgtMapPtr_.reset(new mapDistribute(globalTgtFaces, srcAddress_, cMap));
1077 
1078  if (debug)
1079  {
1080  writeFaceConnectivity(srcPatch, newTgtPatch, srcAddress_);
1081  }
1082  }
1083  else
1084  {
1085  // Calculate AMI interpolation
1086  autoPtr<AMIMethod> AMIPtr
1087  (
1089  (
1090  methodName_,
1091  srcPatch,
1092  tgtPatch,
1093  srcMagSf_,
1094  tgtMagSf_,
1095  triMode_,
1096  reverseTarget_,
1097  requireMatch_ && (lowWeightCorrection_ < 0)
1098  )
1099  );
1100 
1101  AMIPtr->calculate
1102  (
1103  srcAddress_,
1104  srcWeights_,
1105  tgtAddress_,
1106  tgtWeights_
1107  );
1108  }
1109 
1110  // Weight summation and normalisation
1111  sumWeights(*this);
1112  if (report)
1113  {
1114  reportSumWeights(*this);
1115  }
1116  if (requireMatch_)
1117  {
1118  normaliseWeights(*this);
1119  }
1120 
1121  if (debug)
1122  {
1123  Info<< "AMIIPatchToPatchnterpolation :"
1124  << "Constructed addressing and weights" << nl
1125  << " triMode :"
1127  << " singlePatchProc:" << singlePatchProc_ << nl
1128  << " srcMagSf :" << gSum(srcMagSf_) << nl
1129  << " tgtMagSf :" << gSum(tgtMagSf_) << nl
1130  << endl;
1131  }
1132 }
1133 
1134 
1135 void Foam::AMIInterpolation::sumWeights(AMIInterpolation& AMI)
1136 {
1137  sumWeights(AMI.srcWeights_, AMI.srcWeightsSum_);
1138  sumWeights(AMI.tgtWeights_, AMI.tgtWeightsSum_);
1139 }
1140 
1141 
1142 void Foam::AMIInterpolation::sumWeights(PtrList<AMIInterpolation>& AMIs)
1143 {
1144  UPtrList<scalarListList> srcWeights(AMIs.size());
1145  forAll(AMIs, i)
1146  {
1147  srcWeights.set(i, &AMIs[i].srcWeights_);
1148  }
1149 
1150  sumWeights(srcWeights, AMIs[0].srcWeightsSum_);
1151 
1152  for (label i = 1; i < AMIs.size(); ++ i)
1153  {
1154  AMIs[i].srcWeightsSum_ = AMIs[0].srcWeightsSum_;
1155  }
1156 
1157  UPtrList<scalarListList> tgtWeights(AMIs.size());
1158  forAll(AMIs, i)
1159  {
1160  tgtWeights.set(i, &AMIs[i].tgtWeights_);
1161  }
1162 
1163  sumWeights(tgtWeights, AMIs[0].tgtWeightsSum_);
1164 
1165  for (label i = 1; i < AMIs.size(); ++ i)
1166  {
1167  AMIs[i].tgtWeightsSum_ = AMIs[0].tgtWeightsSum_;
1168  }
1169 }
1170 
1171 
1172 void Foam::AMIInterpolation::reportSumWeights(AMIInterpolation& AMI)
1173 {
1174  reportSumWeights
1175  (
1176  AMI.srcMagSf_,
1177  "source",
1178  AMI.srcWeightsSum_,
1179  AMI.lowWeightCorrection_
1180  );
1181 
1182  reportSumWeights
1183  (
1184  AMI.tgtMagSf_,
1185  "target",
1186  AMI.tgtWeightsSum_,
1187  AMI.lowWeightCorrection_
1188  );
1189 }
1190 
1191 
1192 void Foam::AMIInterpolation::normaliseWeights(AMIInterpolation& AMI)
1193 {
1194  normaliseWeights(AMI.srcWeights_, AMI.srcWeightsSum_);
1195  normaliseWeights(AMI.tgtWeights_, AMI.tgtWeightsSum_);
1196 }
1197 
1198 
1199 void Foam::AMIInterpolation::normaliseWeights(UPtrList<AMIInterpolation>& AMIs)
1200 {
1201  UPtrList<scalarListList> srcWeights(AMIs.size());
1202  forAll(AMIs, i)
1203  {
1204  srcWeights.set(i, &AMIs[i].srcWeights_);
1205  }
1206 
1207  normaliseWeights(srcWeights, AMIs[0].srcWeightsSum_);
1208 
1209  UPtrList<scalarListList> tgtWeights(AMIs.size());
1210  forAll(AMIs, i)
1211  {
1212  tgtWeights.set(i, &AMIs[i].tgtWeights_);
1213  }
1214 
1215  normaliseWeights(tgtWeights, AMIs[0].tgtWeightsSum_);
1216 }
1217 
1218 
1221  const primitivePatch& srcPatch,
1222  const primitivePatch& tgtPatch,
1223  const vector& n,
1224  const label tgtFacei,
1225  point& tgtPoint
1226 )
1227 const
1228 {
1229  const pointField& srcPoints = srcPatch.points();
1230 
1231  // Source face addresses that intersect target face tgtFacei
1232  const labelList& addr = tgtAddress_[tgtFacei];
1233 
1234  pointHit nearest;
1235  label nearestFacei = -1;
1236 
1237  forAll(addr, i)
1238  {
1239  const label srcFacei = addr[i];
1240  const face& f = srcPatch[srcFacei];
1241 
1242  const pointHit ray =
1243  f.ray(tgtPoint, n, srcPoints, intersection::algorithm::visible);
1244 
1245  if (ray.hit())
1246  {
1247  tgtPoint = ray.rawPoint();
1248  return srcFacei;
1249  }
1250 
1251  const pointHit near = f.nearestPoint(tgtPoint, srcPoints);
1252 
1253  if (near.distance() < nearest.distance())
1254  {
1255  nearest = near;
1256  nearestFacei = srcFacei;
1257  }
1258  }
1259 
1260  if (nearest.hit() || nearest.eligibleMiss())
1261  {
1262  tgtPoint = nearest.rawPoint();
1263  return nearestFacei;
1264  }
1265 
1266  return -1;
1267 }
1268 
1269 
1272  const primitivePatch& srcPatch,
1273  const primitivePatch& tgtPatch,
1274  const vector& n,
1275  const label srcFacei,
1276  point& srcPoint
1277 )
1278 const
1279 {
1280  const pointField& tgtPoints = tgtPatch.points();
1281 
1282  pointHit nearest;
1283  label nearestFacei = -1;
1284 
1285  // Target face addresses that intersect source face srcFacei
1286  const labelList& addr = srcAddress_[srcFacei];
1287 
1288  forAll(addr, i)
1289  {
1290  const label tgtFacei = addr[i];
1291  const face& f = tgtPatch[tgtFacei];
1292 
1293  const pointHit ray =
1294  f.ray(srcPoint, n, tgtPoints, intersection::algorithm::visible);
1295 
1296  if (ray.hit())
1297  {
1298  srcPoint = ray.rawPoint();
1299  return tgtFacei;
1300  }
1301 
1302  const pointHit near = f.nearestPoint(srcPoint, tgtPoints);
1303 
1304  if (near.distance() < nearest.distance())
1305  {
1306  nearest = near;
1307  nearestFacei = tgtFacei;
1308  }
1309  }
1310 
1311  if (nearest.hit() || nearest.eligibleMiss())
1312  {
1313  srcPoint = nearest.rawPoint();
1314  return nearestFacei;
1315  }
1316 
1317  return -1;
1318 }
1319 
1320 
1323  const primitivePatch& srcPatch,
1324  const primitivePatch& tgtPatch,
1325  const labelListList& srcAddress
1326 )
1327 const
1328 {
1329  OFstream os("faceConnectivity" + Foam::name(Pstream::myProcNo()) + ".obj");
1330 
1331  label ptI = 1;
1332 
1333  forAll(srcAddress, i)
1334  {
1335  const labelList& addr = srcAddress[i];
1336  const point& srcPt = srcPatch.faceCentres()[i];
1337 
1338  forAll(addr, j)
1339  {
1340  label tgtPtI = addr[j];
1341  const point& tgtPt = tgtPatch.faceCentres()[tgtPtI];
1342 
1343  meshTools::writeOBJ(os, srcPt);
1344  meshTools::writeOBJ(os, tgtPt);
1345 
1346  os << "l " << ptI << " " << ptI + 1 << endl;
1347 
1348  ptI += 2;
1349  }
1350  }
1351 }
1352 
1353 
1354 // ************************************************************************* //
static word interpolationMethodToWord(const interpolationMethod &method)
Convert interpolationMethod to word representation.
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
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const Field< PointType > & faceCentres() const
Return face centres for patch.
Type gMin(const FieldField< Field, Type > &f)
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
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:174
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.
const labelListList & subMap() const
From subsetted data back to original data.
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
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 Field< PointType > & localPoints() const
Return pointField of points in patch.
static const NamedEnum< triangulationMode, 2 > triangulationModeNames_
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
interpolationMethod
Enumeration specifying interpolation method.
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
AMIInterpolation(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const faceAreaIntersect::triangulationMode &triMode, const bool requireMatch=true, const interpolationMethod &method=imFaceAreaWeight, const scalar lowWeightCorrection=-1, const bool reverseTarget=false, const bool report=true)
Construct from components.
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
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.
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
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
label constructSize() const
Constructed data size.
void setSize(const label)
Reset size of List.
Definition: List.C:281
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
Class containing processor-to-processor mapping information.
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
Input from memory buffer stream.
Definition: IStringStream.H:49
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
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...
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
messageStream Info
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Helper class for list to append y onto the end of x.
Definition: ListOps.H:259
static interpolationMethod wordTointerpolationMethod(const word &method)
Convert word to interpolationMethod.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
static void triangulate(const face &f, const pointField &points, const triangulationMode &triMode, faceList &faceTris)
Triangulate a face using the given triangulation mode.
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.
Namespace for OpenFOAM.