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-2018 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 template<class SourcePatch, class TargetPatch>
37 (
38  const interpolationMethod& im
39 )
40 {
41  word method = "unknown-interpolationMethod";
42 
43  switch (im)
44  {
45  case imDirect:
46  {
47  method = "directAMI";
48  break;
49  }
50  case imMapNearest:
51  {
52  method = "mapNearestAMI";
53  break;
54  }
55  case imFaceAreaWeight:
56  {
57  method = "faceAreaWeightAMI";
58  break;
59  }
60  case imPartialFaceAreaWeight:
61  {
62  method = "partialFaceAreaWeightAMI";
63  break;
64  }
65  case imSweptFaceAreaWeight:
66  {
67  method = "sweptFaceAreaWeightAMI";
68  break;
69  }
70  default:
71  {
73  << "Unhandled interpolationMethod enumeration " << method
74  << abort(FatalError);
75  }
76  }
77 
78  return method;
79 }
80 
81 
82 template<class SourcePatch, class TargetPatch>
85 (
86  const word& im
87 )
88 {
89  interpolationMethod method = imDirect;
90 
91  wordList methods
92  (
94  (
95  "("
96  "directAMI "
97  "mapNearestAMI "
98  "faceAreaWeightAMI "
99  "partialFaceAreaWeightAMI "
100  "sweptFaceAreaWeightAMI"
101  ")"
102  )()
103  );
104 
105  if (im == "directAMI")
106  {
107  method = imDirect;
108  }
109  else if (im == "mapNearestAMI")
110  {
111  method = imMapNearest;
112  }
113  else if (im == "faceAreaWeightAMI")
114  {
115  method = imFaceAreaWeight;
116  }
117  else if (im == "partialFaceAreaWeightAMI")
118  {
119  method = imPartialFaceAreaWeight;
120  }
121  else if (im == "sweptFaceAreaWeightAMI")
122  {
123  method = imSweptFaceAreaWeight;
124  }
125  else
126  {
128  << "Invalid interpolationMethod " << im
129  << ". Valid methods are:" << methods
130  << exit(FatalError);
131  }
132 
133  return method;
134 }
135 
136 
137 template<class SourcePatch, class TargetPatch>
138 template<class Patch>
141 (
142  const Patch& patch,
144 )
145 {
146  tmp<scalarField> tResult(new scalarField(patch.size(), Zero));
147  scalarField& result = tResult.ref();
148 
149  const pointField& patchPoints = patch.localPoints();
150 
151  faceList patchFaceTris;
152 
153  forAll(result, patchFacei)
154  {
155  faceAreaIntersect::triangulate
156  (
157  patch.localFaces()[patchFacei],
158  patchPoints,
159  triMode,
160  patchFaceTris
161  );
162 
163  forAll(patchFaceTris, i)
164  {
165  result[patchFacei] +=
167  (
168  patchPoints[patchFaceTris[i][0]],
169  patchPoints[patchFaceTris[i][1]],
170  patchPoints[patchFaceTris[i][2]]
171  ).mag();
172  }
173  }
174 
175  return tResult;
176 }
177 
178 
179 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
180 
181 template<class SourcePatch, class TargetPatch>
183 (
184  const searchableSurface& surf,
185  pointField& pts
186 ) const
187 {
188  if (debug)
189  {
190  Info<< "AMI: projecting points to surface" << endl;
191  }
192 
194 
195  surf.findNearest(pts, scalarField(pts.size(), great), nearInfo);
196 
197  label nMiss = 0;
198  forAll(nearInfo, i)
199  {
200  const pointIndexHit& pi = nearInfo[i];
201 
202  if (pi.hit())
203  {
204  pts[i] = pi.hitPoint();
205  }
206  else
207  {
208  pts[i] = pts[i];
209  nMiss++;
210  }
211  }
212 
213  if (nMiss > 0)
214  {
216  << "Error projecting points to surface: "
217  << nMiss << " faces could not be determined"
218  << abort(FatalError);
219  }
220 }
221 
222 
223 template<class SourcePatch, class TargetPatch>
225 (
226  const scalarListList& wght,
227  scalarField& wghtSum
228 )
229 {
230  wghtSum.setSize(wght.size());
231  wghtSum = Zero;
232 
233  forAll(wght, facei)
234  {
235  wghtSum[facei] = sum(wght[facei]);
236  }
237 }
238 
239 
240 template<class SourcePatch, class TargetPatch>
242 (
243  const UPtrList<scalarListList>& wghts,
244  scalarField& wghtSum
245 )
246 {
247  wghtSum.setSize(wghts[0].size());
248  wghtSum = Zero;
249 
250  forAll(wghts[0], facei)
251  {
252  forAll(wghts, wghtsi)
253  {
254  forAll(wghts[wghtsi][facei], i)
255  {
256  wghtSum[facei] += wghts[wghtsi][facei][i];
257  }
258  }
259  }
260 }
261 
262 
263 template<class SourcePatch, class TargetPatch>
265 (
266  const scalarField& patchAreas,
267  const word& patchName,
268  const scalarField& wghtSum,
269  const scalar lowWeightTol
270 )
271 {
272  if (returnReduce(wghtSum.size(), sumOp<label>()) == 0)
273  {
274  return;
275  }
276 
277  label nLowWeight = 0;
278  forAll(wghtSum, facei)
279  {
280  if (wghtSum[facei] < lowWeightTol)
281  {
282  ++ nLowWeight;
283  }
284  }
285  reduce(nLowWeight, sumOp<label>());
286 
287  Info<< indent << "AMI: Patch " << patchName
288  << " sum(weights) min/max/average = " << gMin(wghtSum) << ", "
289  << gMax(wghtSum) << ", "
290  << gSum(wghtSum*patchAreas)/gSum(patchAreas) << endl;
291 
292  if (nLowWeight)
293  {
294  Info<< indent << "AMI: Patch " << patchName << " identified "
295  << nLowWeight << " faces with weights less than "
296  << lowWeightTol << endl;
297  }
298 }
299 
300 
301 template<class SourcePatch, class TargetPatch>
303 (
304  scalarListList& wght,
305  const scalarField& wghtSum
306 )
307 {
308  forAll(wghtSum, facei)
309  {
310  scalarList& w = wght[facei];
311 
312  forAll(w, i)
313  {
314  w[i] /= wghtSum[facei];
315  }
316  }
317 }
318 
319 
320 template<class SourcePatch, class TargetPatch>
322 (
324  const scalarField& wghtSum
325 )
326 {
327  forAll(wghtSum, facei)
328  {
329  forAll(wghts, wghtsi)
330  {
331  scalarList& w = wghts[wghtsi][facei];
332 
333  forAll(w, i)
334  {
335  w[i] /= wghtSum[facei];
336  }
337  }
338  }
339 }
340 
341 
342 template<class SourcePatch, class TargetPatch>
344 (
345  const autoPtr<mapDistribute>& targetMapPtr,
346  const scalarField& fineSrcMagSf,
347  const labelListList& fineSrcAddress,
348  const scalarListList& fineSrcWeights,
349 
350  const labelList& sourceRestrictAddressing,
351  const labelList& targetRestrictAddressing,
352 
353  scalarField& srcMagSf,
354  labelListList& srcAddress,
355  scalarListList& srcWeights,
356  scalarField& srcWeightsSum,
357  autoPtr<mapDistribute>& tgtMap
358 )
359 {
360  label sourceCoarseSize =
361  (
362  sourceRestrictAddressing.size()
363  ? max(sourceRestrictAddressing)+1
364  : 0
365  );
366 
367  label targetCoarseSize =
368  (
369  targetRestrictAddressing.size()
370  ? max(targetRestrictAddressing)+1
371  : 0
372  );
373 
374  // Agglomerate face areas
375  {
376  srcMagSf.setSize(sourceRestrictAddressing.size(), 0.0);
377  forAll(sourceRestrictAddressing, facei)
378  {
379  label coarseFacei = sourceRestrictAddressing[facei];
380  srcMagSf[coarseFacei] += fineSrcMagSf[facei];
381  }
382  }
383 
384 
385  // Agglomerate weights and indices
386  if (targetMapPtr.valid())
387  {
388  const mapDistribute& map = targetMapPtr();
389 
390  // Get all restriction addressing.
391  labelList allRestrict(targetRestrictAddressing);
392  map.distribute(allRestrict);
393 
394  // So now we have agglomeration of the target side in
395  // allRestrict:
396  // 0..size-1 : local agglomeration (= targetRestrictAddressing)
397  // size.. : agglomeration data from other processors
398 
399  labelListList tgtSubMap(Pstream::nProcs());
400 
401  // Local subMap is just identity
402  {
403  tgtSubMap[Pstream::myProcNo()] = identity(targetCoarseSize);
404  }
405 
406  forAll(map.subMap(), proci)
407  {
408  if (proci != Pstream::myProcNo())
409  {
410  // Combine entries that point to the same coarse element. All
411  // the elements refer to local data so index into
412  // targetRestrictAddressing or allRestrict (since the same
413  // for local data).
414  const labelList& elems = map.subMap()[proci];
415  labelList& newSubMap = tgtSubMap[proci];
416  newSubMap.setSize(elems.size());
417 
418  labelList oldToNew(targetCoarseSize, -1);
419  label newI = 0;
420 
421  forAll(elems, i)
422  {
423  label fineElem = elems[i];
424  label coarseElem = allRestrict[fineElem];
425  if (oldToNew[coarseElem] == -1)
426  {
427  oldToNew[coarseElem] = newI;
428  newSubMap[newI] = coarseElem;
429  newI++;
430  }
431  }
432  newSubMap.setSize(newI);
433  }
434  }
435 
436  // Reconstruct constructMap by combining entries. Note that order
437  // of handing out indices should be the same as loop above to compact
438  // the sending map
439 
440  labelListList tgtConstructMap(Pstream::nProcs());
441  labelList tgtCompactMap;
442 
443  // Local constructMap is just identity
444  {
445  tgtConstructMap[Pstream::myProcNo()] =
446  identity(targetCoarseSize);
447 
448  tgtCompactMap = targetRestrictAddressing;
449  }
450  tgtCompactMap.setSize(map.constructSize());
451  label compactI = targetCoarseSize;
452 
453  // Compact data from other processors
454  forAll(map.constructMap(), proci)
455  {
456  if (proci != Pstream::myProcNo())
457  {
458  // Combine entries that point to the same coarse element. All
459  // elements now are remote data so we cannot use any local
460  // data here - use allRestrict instead.
461  const labelList& elems = map.constructMap()[proci];
462 
463  labelList& newConstructMap = tgtConstructMap[proci];
464  newConstructMap.setSize(elems.size());
465 
466  if (elems.size())
467  {
468  // Get the maximum target coarse size for this set of
469  // received data.
470  label remoteTargetCoarseSize = labelMin;
471  forAll(elems, i)
472  {
473  remoteTargetCoarseSize = max
474  (
475  remoteTargetCoarseSize,
476  allRestrict[elems[i]]
477  );
478  }
479  remoteTargetCoarseSize += 1;
480 
481  // Combine locally data coming from proci
482  labelList oldToNew(remoteTargetCoarseSize, -1);
483  label newI = 0;
484 
485  forAll(elems, i)
486  {
487  label fineElem = elems[i];
488  // fineElem now points to section from proci
489  label coarseElem = allRestrict[fineElem];
490  if (oldToNew[coarseElem] == -1)
491  {
492  oldToNew[coarseElem] = newI;
493  tgtCompactMap[fineElem] = compactI;
494  newConstructMap[newI] = compactI++;
495  newI++;
496  }
497  else
498  {
499  // Get compact index
500  label compactI = oldToNew[coarseElem];
501  tgtCompactMap[fineElem] = newConstructMap[compactI];
502  }
503  }
504  newConstructMap.setSize(newI);
505  }
506  }
507  }
508 
509 
510  srcAddress.setSize(sourceCoarseSize);
511  srcWeights.setSize(sourceCoarseSize);
512 
513  forAll(fineSrcAddress, facei)
514  {
515  // All the elements contributing to facei. Are slots in
516  // mapDistribute'd data.
517  const labelList& elems = fineSrcAddress[facei];
518  const scalarList& weights = fineSrcWeights[facei];
519  const scalar fineArea = fineSrcMagSf[facei];
520 
521  label coarseFacei = sourceRestrictAddressing[facei];
522 
523  const scalar coarseArea = srcMagSf[coarseFacei];
524 
525  labelList& newElems = srcAddress[coarseFacei];
526  scalarList& newWeights = srcWeights[coarseFacei];
527 
528  forAll(elems, i)
529  {
530  label elemI = elems[i];
531  label coarseElemI = tgtCompactMap[elemI];
532 
533  label index = findIndex(newElems, coarseElemI);
534  if (index == -1)
535  {
536  newElems.append(coarseElemI);
537  newWeights.append(fineArea/coarseArea*weights[i]);
538  }
539  else
540  {
541  newWeights[index] += fineArea/coarseArea*weights[i];
542  }
543  }
544  }
545 
546  tgtMap.reset
547  (
548  new mapDistribute
549  (
550  compactI,
551  tgtSubMap.xfer(),
552  tgtConstructMap.xfer()
553  )
554  );
555  }
556  else
557  {
558  srcAddress.setSize(sourceCoarseSize);
559  srcWeights.setSize(sourceCoarseSize);
560 
561  forAll(fineSrcAddress, facei)
562  {
563  // All the elements contributing to facei. Are slots in
564  // mapDistribute'd data.
565  const labelList& elems = fineSrcAddress[facei];
566  const scalarList& weights = fineSrcWeights[facei];
567  const scalar fineArea = fineSrcMagSf[facei];
568 
569  label coarseFacei = sourceRestrictAddressing[facei];
570 
571  const scalar coarseArea = srcMagSf[coarseFacei];
572 
573  labelList& newElems = srcAddress[coarseFacei];
574  scalarList& newWeights = srcWeights[coarseFacei];
575 
576  forAll(elems, i)
577  {
578  label elemI = elems[i];
579  label coarseElemI = targetRestrictAddressing[elemI];
580 
581  label index = findIndex(newElems, coarseElemI);
582  if (index == -1)
583  {
584  newElems.append(coarseElemI);
585  newWeights.append(fineArea/coarseArea*weights[i]);
586  }
587  else
588  {
589  newWeights[index] += fineArea/coarseArea*weights[i];
590  }
591  }
592  }
593  }
594 }
595 
596 
597 template<class SourcePatch, class TargetPatch>
599 (
600  const SourcePatch& srcPatch,
601  const TargetPatch& tgtPatch,
602  const autoPtr<searchableSurface>& surfPtr,
603  const bool report
604 )
605 {
606  if (surfPtr.valid())
607  {
608  // Create new patches for source and target
609  pointField srcPoints = srcPatch.points();
610  SourcePatch srcPatch0
611  (
613  (
614  srcPatch,
615  srcPatch.size(),
616  0
617  ),
618  srcPoints
619  );
620 
621  if (debug)
622  {
623  OFstream os("amiSrcPoints.obj");
624  forAll(srcPoints, i)
625  {
626  meshTools::writeOBJ(os, srcPoints[i]);
627  }
628  }
629 
630  pointField tgtPoints = tgtPatch.points();
631  TargetPatch tgtPatch0
632  (
634  (
635  tgtPatch,
636  tgtPatch.size(),
637  0
638  ),
639  tgtPoints
640  );
641 
642  if (debug)
643  {
644  OFstream os("amiTgtPoints.obj");
645  forAll(tgtPoints, i)
646  {
647  meshTools::writeOBJ(os, tgtPoints[i]);
648  }
649  }
650 
651 
652  // Map source and target patches onto projection surface
653  projectPointsToSurface(surfPtr(), srcPoints);
654  projectPointsToSurface(surfPtr(), tgtPoints);
655 
656 
657  // Calculate AMI interpolation
658  update(srcPatch0, tgtPatch0, report);
659  }
660  else
661  {
662  update(srcPatch, tgtPatch, report);
663  }
664 }
665 
666 
667 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
668 
669 template<class SourcePatch, class TargetPatch>
671 (
672  const SourcePatch& srcPatch,
673  const TargetPatch& tgtPatch,
675  const bool requireMatch,
676  const interpolationMethod& method,
677  const scalar lowWeightCorrection,
678  const bool reverseTarget,
679  const bool report
680 )
681 :
682  methodName_(interpolationMethodToWord(method)),
683  reverseTarget_(reverseTarget),
684  requireMatch_(requireMatch),
685  singlePatchProc_(-999),
686  lowWeightCorrection_(lowWeightCorrection),
687  srcAddress_(),
688  srcWeights_(),
689  srcWeightsSum_(),
690  tgtAddress_(),
691  tgtWeights_(),
692  tgtWeightsSum_(),
693  triMode_(triMode),
694  srcMapPtr_(nullptr),
695  tgtMapPtr_(nullptr)
696 {
697  update(srcPatch, tgtPatch, report);
698 }
699 
700 
701 template<class SourcePatch, class TargetPatch>
703 (
704  const SourcePatch& srcPatch,
705  const TargetPatch& tgtPatch,
707  const bool requireMatch,
708  const word& methodName,
709  const scalar lowWeightCorrection,
710  const bool reverseTarget,
711  const bool report
712 )
713 :
714  methodName_(methodName),
715  reverseTarget_(reverseTarget),
716  requireMatch_(requireMatch),
717  singlePatchProc_(-999),
718  lowWeightCorrection_(lowWeightCorrection),
719  srcAddress_(),
720  srcWeights_(),
721  srcWeightsSum_(),
722  tgtAddress_(),
723  tgtWeights_(),
724  tgtWeightsSum_(),
725  triMode_(triMode),
726  srcMapPtr_(nullptr),
727  tgtMapPtr_(nullptr)
728 {
729  update(srcPatch, tgtPatch, report);
730 }
731 
732 
733 template<class SourcePatch, class TargetPatch>
735 (
736  const SourcePatch& srcPatch,
737  const TargetPatch& tgtPatch,
738  const autoPtr<searchableSurface>& surfPtr,
740  const bool requireMatch,
741  const interpolationMethod& method,
742  const scalar lowWeightCorrection,
743  const bool reverseTarget,
744  const bool report
745 )
746 :
747  methodName_(interpolationMethodToWord(method)),
748  reverseTarget_(reverseTarget),
749  requireMatch_(requireMatch),
750  singlePatchProc_(-999),
751  lowWeightCorrection_(lowWeightCorrection),
752  srcAddress_(),
753  srcWeights_(),
754  srcWeightsSum_(),
755  tgtAddress_(),
756  tgtWeights_(),
757  tgtWeightsSum_(),
758  triMode_(triMode),
759  srcMapPtr_(nullptr),
760  tgtMapPtr_(nullptr)
761 {
762  constructFromSurface(srcPatch, tgtPatch, surfPtr, report);
763 }
764 
765 
766 template<class SourcePatch, class TargetPatch>
768 (
769  const SourcePatch& srcPatch,
770  const TargetPatch& tgtPatch,
771  const autoPtr<searchableSurface>& surfPtr,
773  const bool requireMatch,
774  const word& methodName,
775  const scalar lowWeightCorrection,
776  const bool reverseTarget,
777  const bool report
778 )
779 :
780  methodName_(methodName),
781  reverseTarget_(reverseTarget),
782  requireMatch_(requireMatch),
783  singlePatchProc_(-999),
784  lowWeightCorrection_(lowWeightCorrection),
785  srcAddress_(),
786  srcWeights_(),
787  srcWeightsSum_(),
788  tgtAddress_(),
789  tgtWeights_(),
790  tgtWeightsSum_(),
791  triMode_(triMode),
792  srcMapPtr_(nullptr),
793  tgtMapPtr_(nullptr)
794 {
795  constructFromSurface(srcPatch, tgtPatch, surfPtr, report);
796 }
797 
798 
799 template<class SourcePatch, class TargetPatch>
801 (
803  const labelList& sourceRestrictAddressing,
804  const labelList& targetRestrictAddressing,
805  const bool report
806 )
807 :
808  methodName_(fineAMI.methodName_),
809  reverseTarget_(fineAMI.reverseTarget_),
810  requireMatch_(fineAMI.requireMatch_),
811  singlePatchProc_(fineAMI.singlePatchProc_),
812  lowWeightCorrection_(-1.0),
813  srcAddress_(),
814  srcWeights_(),
815  srcWeightsSum_(),
816  tgtAddress_(),
817  tgtWeights_(),
818  tgtWeightsSum_(),
819  triMode_(fineAMI.triMode_),
820  srcMapPtr_(nullptr),
821  tgtMapPtr_(nullptr)
822 {
823  label sourceCoarseSize =
824  (
825  sourceRestrictAddressing.size()
826  ? max(sourceRestrictAddressing)+1
827  : 0
828  );
829 
830  label neighbourCoarseSize =
831  (
832  targetRestrictAddressing.size()
833  ? max(targetRestrictAddressing)+1
834  : 0
835  );
836 
837  if (debug & 2)
838  {
839  Pout<< "AMI: Creating addressing and weights as agglomeration of AMI :"
840  << " source:" << fineAMI.srcAddress().size()
841  << " target:" << fineAMI.tgtAddress().size()
842  << " coarse source size:" << sourceCoarseSize
843  << " neighbour source size:" << neighbourCoarseSize
844  << endl;
845  }
846 
847  if
848  (
849  fineAMI.srcAddress().size() != sourceRestrictAddressing.size()
850  || fineAMI.tgtAddress().size() != targetRestrictAddressing.size()
851  )
852  {
854  << "Size mismatch." << nl
855  << "Source patch size:" << fineAMI.srcAddress().size() << nl
856  << "Source agglomeration size:"
857  << sourceRestrictAddressing.size() << nl
858  << "Target patch size:" << fineAMI.tgtAddress().size() << nl
859  << "Target agglomeration size:"
860  << targetRestrictAddressing.size()
861  << exit(FatalError);
862  }
863 
864  // Agglomerate addresses and weights
865  agglomerate
866  (
867  fineAMI.tgtMapPtr_,
868  fineAMI.srcMagSf(),
869  fineAMI.srcAddress(),
870  fineAMI.srcWeights(),
871 
872  sourceRestrictAddressing,
873  targetRestrictAddressing,
874 
875  srcMagSf_,
876  srcAddress_,
877  srcWeights_,
878  srcWeightsSum_,
879  tgtMapPtr_
880  );
881 
882  agglomerate
883  (
884  fineAMI.srcMapPtr_,
885  fineAMI.tgtMagSf(),
886  fineAMI.tgtAddress(),
887  fineAMI.tgtWeights(),
888 
889  targetRestrictAddressing,
890  sourceRestrictAddressing,
891 
892  tgtMagSf_,
893  tgtAddress_,
894  tgtWeights_,
895  tgtWeightsSum_,
896  srcMapPtr_
897  );
898 
899  // Weight summation and normalisation
900  sumWeights(*this);
901  if (report)
902  {
903  reportSumWeights(*this);
904  }
905  if (requireMatch_)
906  {
907  normaliseWeights(*this);
908  }
909 }
910 
911 
912 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
913 
914 template<class SourcePatch, class TargetPatch>
916 {}
917 
918 
919 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
920 
921 template<class SourcePatch, class TargetPatch>
923 (
924  const SourcePatch& srcPatch,
925  const TargetPatch& tgtPatch,
926  const bool report
927 )
928 {
929  label srcTotalSize = returnReduce(srcPatch.size(), sumOp<label>());
930  label tgtTotalSize = returnReduce(tgtPatch.size(), sumOp<label>());
931 
932  if (srcTotalSize == 0)
933  {
934  if (debug)
935  {
936  Info<< "AMI: no source faces present - no addressing constructed"
937  << endl;
938  }
939 
940  return;
941  }
942 
943  if (report)
944  {
945  Info<< indent
946  << "AMI: Creating addressing and weights between "
947  << srcTotalSize << " source faces and "
948  << tgtTotalSize << " target faces"
949  << endl;
950  }
951 
952  // Calculate face areas
953  srcMagSf_ = patchMagSf(srcPatch, triMode_);
954  tgtMagSf_ = patchMagSf(tgtPatch, triMode_);
955 
956  // Calculate if patches present on multiple processors
957  singlePatchProc_ = calcDistribution(srcPatch, tgtPatch);
958 
959  if (singlePatchProc_ == -1)
960  {
961  // Convert local addressing to global addressing
962  globalIndex globalSrcFaces(srcPatch.size());
963  globalIndex globalTgtFaces(tgtPatch.size());
964 
965  // Create processor map of overlapping faces. This map gets
966  // (possibly remote) faces from the tgtPatch such that they (together)
967  // cover all of the srcPatch
968  autoPtr<mapDistribute> mapPtr = calcProcMap(srcPatch, tgtPatch);
969  const mapDistribute& map = mapPtr();
970 
971  // Create new target patch that fully encompasses source patch
972 
973  // Faces and points
974  faceList newTgtFaces;
975  pointField newTgtPoints;
976 
977  // Original faces from tgtPatch (in globalIndexing since might be
978  // remote)
979  labelList tgtFaceIDs;
980  distributeAndMergePatches
981  (
982  map,
983  tgtPatch,
984  globalTgtFaces,
985  newTgtFaces,
986  newTgtPoints,
987  tgtFaceIDs
988  );
989 
990  TargetPatch
991  newTgtPatch
992  (
994  (
995  newTgtFaces,
996  newTgtFaces.size()
997  ),
998  newTgtPoints
999  );
1000  scalarField newTgtMagSf(patchMagSf(newTgtPatch, triMode_));
1001 
1002  // Calculate AMI interpolation
1004  (
1006  (
1007  methodName_,
1008  srcPatch,
1009  newTgtPatch,
1010  srcMagSf_,
1011  newTgtMagSf,
1012  triMode_,
1013  reverseTarget_,
1014  requireMatch_ && (lowWeightCorrection_ < 0)
1015  )
1016  );
1017 
1018  AMIPtr->calculate
1019  (
1020  srcAddress_,
1021  srcWeights_,
1022  tgtAddress_,
1023  tgtWeights_
1024  );
1025 
1026  // Now
1027  // ~~~
1028  // srcAddress_ : per srcPatch face a list of the newTgtPatch (not
1029  // tgtPatch) faces it overlaps
1030  // tgtAddress_ : per newTgtPatch (not tgtPatch) face a list of the
1031  // srcPatch faces it overlaps
1032 
1033 
1034  // Rework newTgtPatch indices into globalIndices of tgtPatch
1035  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1036 
1037 
1038  forAll(srcAddress_, i)
1039  {
1040  labelList& addressing = srcAddress_[i];
1041  forAll(addressing, addrI)
1042  {
1043  addressing[addrI] = tgtFaceIDs[addressing[addrI]];
1044  }
1045  }
1046 
1047  forAll(tgtAddress_, i)
1048  {
1049  labelList& addressing = tgtAddress_[i];
1050  forAll(addressing, addrI)
1051  {
1052  addressing[addrI] = globalSrcFaces.toGlobal(addressing[addrI]);
1053  }
1054  }
1055 
1056  // Send data back to originating procs. Note that contributions
1057  // from different processors get added (ListAppendEqOp)
1058 
1059  mapDistributeBase::distribute
1060  (
1061  Pstream::commsTypes::nonBlocking,
1062  List<labelPair>(),
1063  tgtPatch.size(),
1064  map.constructMap(),
1065  false, // has flip
1066  map.subMap(),
1067  false, // has flip
1068  tgtAddress_,
1070  flipOp(), // flip operation
1071  labelList()
1072  );
1073 
1074  mapDistributeBase::distribute
1075  (
1076  Pstream::commsTypes::nonBlocking,
1077  List<labelPair>(),
1078  tgtPatch.size(),
1079  map.constructMap(),
1080  false,
1081  map.subMap(),
1082  false,
1083  tgtWeights_,
1085  flipOp(),
1086  scalarList()
1087  );
1088 
1089  // Cache maps and reset addresses
1090  List<Map<label>> cMap;
1091  srcMapPtr_.reset(new mapDistribute(globalSrcFaces, tgtAddress_, cMap));
1092  tgtMapPtr_.reset(new mapDistribute(globalTgtFaces, srcAddress_, cMap));
1093 
1094  if (debug)
1095  {
1096  writeFaceConnectivity(srcPatch, newTgtPatch, srcAddress_);
1097  }
1098  }
1099  else
1100  {
1101  // Calculate AMI interpolation
1103  (
1105  (
1106  methodName_,
1107  srcPatch,
1108  tgtPatch,
1109  srcMagSf_,
1110  tgtMagSf_,
1111  triMode_,
1112  reverseTarget_,
1113  requireMatch_ && (lowWeightCorrection_ < 0)
1114  )
1115  );
1116 
1117  AMIPtr->calculate
1118  (
1119  srcAddress_,
1120  srcWeights_,
1121  tgtAddress_,
1122  tgtWeights_
1123  );
1124  }
1125 
1126  // Weight summation and normalisation
1127  sumWeights(*this);
1128  if (report)
1129  {
1130  reportSumWeights(*this);
1131  }
1132  if (requireMatch_)
1133  {
1134  normaliseWeights(*this);
1135  }
1136 
1137  if (debug)
1138  {
1139  Info<< "AMIInterpolation : Constructed addressing and weights" << nl
1140  << " triMode :"
1141  << faceAreaIntersect::triangulationModeNames_[triMode_] << nl
1142  << " singlePatchProc:" << singlePatchProc_ << nl
1143  << " srcMagSf :" << gSum(srcMagSf_) << nl
1144  << " tgtMagSf :" << gSum(tgtMagSf_) << nl
1145  << endl;
1146  }
1147 }
1148 
1149 
1150 template<class SourcePatch, class TargetPatch>
1154 )
1155 {
1156  sumWeights(AMI.srcWeights_, AMI.srcWeightsSum_);
1157  sumWeights(AMI.tgtWeights_, AMI.tgtWeightsSum_);
1158 }
1159 
1160 
1161 template<class SourcePatch, class TargetPatch>
1165 )
1166 {
1167  UPtrList<scalarListList> srcWeights(AMIs.size());
1168  forAll(AMIs, i)
1169  {
1170  srcWeights.set(i, &AMIs[i].srcWeights_);
1171  }
1172 
1173  sumWeights(srcWeights, AMIs[0].srcWeightsSum_);
1174 
1175  for (label i = 1; i < AMIs.size(); ++ i)
1176  {
1177  AMIs[i].srcWeightsSum_ = AMIs[0].srcWeightsSum_;
1178  }
1179 
1180  UPtrList<scalarListList> tgtWeights(AMIs.size());
1181  forAll(AMIs, i)
1182  {
1183  tgtWeights.set(i, &AMIs[i].tgtWeights_);
1184  }
1185 
1186  sumWeights(tgtWeights, AMIs[0].tgtWeightsSum_);
1187 
1188  for (label i = 1; i < AMIs.size(); ++ i)
1189  {
1190  AMIs[i].tgtWeightsSum_ = AMIs[0].tgtWeightsSum_;
1191  }
1192 }
1193 
1194 
1195 template<class SourcePatch, class TargetPatch>
1199 )
1200 {
1201  reportSumWeights
1202  (
1203  AMI.srcMagSf_,
1204  "source",
1205  AMI.srcWeightsSum_,
1206  AMI.lowWeightCorrection_
1207  );
1208 
1209  reportSumWeights
1210  (
1211  AMI.tgtMagSf_,
1212  "target",
1213  AMI.tgtWeightsSum_,
1214  AMI.lowWeightCorrection_
1215  );
1216 }
1217 
1218 
1219 template<class SourcePatch, class TargetPatch>
1223 )
1224 {
1225  normaliseWeights(AMI.srcWeights_, AMI.srcWeightsSum_);
1226  normaliseWeights(AMI.tgtWeights_, AMI.tgtWeightsSum_);
1227 }
1228 
1229 
1230 template<class SourcePatch, class TargetPatch>
1234 )
1235 {
1236  UPtrList<scalarListList> srcWeights(AMIs.size());
1237  forAll(AMIs, i)
1238  {
1239  srcWeights.set(i, &AMIs[i].srcWeights_);
1240  }
1241 
1242  normaliseWeights(srcWeights, AMIs[0].srcWeightsSum_);
1243 
1244  UPtrList<scalarListList> tgtWeights(AMIs.size());
1245  forAll(AMIs, i)
1246  {
1247  tgtWeights.set(i, &AMIs[i].tgtWeights_);
1248  }
1249 
1250  normaliseWeights(tgtWeights, AMIs[0].tgtWeightsSum_);
1251 }
1252 
1253 
1254 template<class SourcePatch, class TargetPatch>
1255 template<class Type, class CombineOp>
1258  const UList<Type>& fld,
1259  const CombineOp& cop,
1260  List<Type>& result,
1261  const UList<Type>& defaultValues
1262 ) const
1263 {
1264  if (fld.size() != srcAddress_.size())
1265  {
1267  << "Supplied field size is not equal to source patch size" << nl
1268  << " source patch = " << srcAddress_.size() << nl
1269  << " target patch = " << tgtAddress_.size() << nl
1270  << " supplied field = " << fld.size()
1271  << abort(FatalError);
1272  }
1273 
1274  if (lowWeightCorrection_ > 0)
1275  {
1276  if (defaultValues.size() != tgtAddress_.size())
1277  {
1279  << "Employing default values when sum of weights falls below "
1280  << lowWeightCorrection_
1281  << " but supplied default field size is not equal to target "
1282  << "patch size" << nl
1283  << " default values = " << defaultValues.size() << nl
1284  << " target patch = " << tgtAddress_.size() << nl
1285  << abort(FatalError);
1286  }
1287  }
1288 
1289  result.setSize(tgtAddress_.size());
1290 
1291  if (singlePatchProc_ == -1)
1292  {
1293  const mapDistribute& map = srcMapPtr_();
1294 
1295  List<Type> work(fld);
1296  map.distribute(work);
1297 
1298  forAll(result, facei)
1299  {
1300  if (tgtWeightsSum_[facei] < lowWeightCorrection_)
1301  {
1302  result[facei] = defaultValues[facei];
1303  }
1304  else
1305  {
1306  const labelList& faces = tgtAddress_[facei];
1307  const scalarList& weights = tgtWeights_[facei];
1308 
1309  forAll(faces, i)
1310  {
1311  cop(result[facei], facei, work[faces[i]], weights[i]);
1312  }
1313  }
1314  }
1315  }
1316  else
1317  {
1318  forAll(result, facei)
1319  {
1320  if (tgtWeightsSum_[facei] < lowWeightCorrection_)
1321  {
1322  result[facei] = defaultValues[facei];
1323  }
1324  else
1325  {
1326  const labelList& faces = tgtAddress_[facei];
1327  const scalarList& weights = tgtWeights_[facei];
1328 
1329  forAll(faces, i)
1330  {
1331  cop(result[facei], facei, fld[faces[i]], weights[i]);
1332  }
1333  }
1334  }
1335  }
1336 }
1337 
1338 
1339 template<class SourcePatch, class TargetPatch>
1340 template<class Type, class CombineOp>
1343  const UList<Type>& fld,
1344  const CombineOp& cop,
1345  List<Type>& result,
1346  const UList<Type>& defaultValues
1347 ) const
1348 {
1349  if (fld.size() != tgtAddress_.size())
1350  {
1352  << "Supplied field size is not equal to target patch size" << nl
1353  << " source patch = " << srcAddress_.size() << nl
1354  << " target patch = " << tgtAddress_.size() << nl
1355  << " supplied field = " << fld.size()
1356  << abort(FatalError);
1357  }
1358 
1359  if (lowWeightCorrection_ > 0)
1360  {
1361  if (defaultValues.size() != srcAddress_.size())
1362  {
1364  << "Employing default values when sum of weights falls below "
1365  << lowWeightCorrection_
1366  << " but supplied default field size is not equal to target "
1367  << "patch size" << nl
1368  << " default values = " << defaultValues.size() << nl
1369  << " source patch = " << srcAddress_.size() << nl
1370  << abort(FatalError);
1371  }
1372  }
1373 
1374  result.setSize(srcAddress_.size());
1375 
1376  if (singlePatchProc_ == -1)
1377  {
1378  const mapDistribute& map = tgtMapPtr_();
1379 
1380  List<Type> work(fld);
1381  map.distribute(work);
1382 
1383  forAll(result, facei)
1384  {
1385  if (srcWeightsSum_[facei] < lowWeightCorrection_)
1386  {
1387  result[facei] = defaultValues[facei];
1388  }
1389  else
1390  {
1391  const labelList& faces = srcAddress_[facei];
1392  const scalarList& weights = srcWeights_[facei];
1393 
1394  forAll(faces, i)
1395  {
1396  cop(result[facei], facei, work[faces[i]], weights[i]);
1397  }
1398  }
1399  }
1400  }
1401  else
1402  {
1403  forAll(result, facei)
1404  {
1405  if (srcWeightsSum_[facei] < lowWeightCorrection_)
1406  {
1407  result[facei] = defaultValues[facei];
1408  }
1409  else
1410  {
1411  const labelList& faces = srcAddress_[facei];
1412  const scalarList& weights = srcWeights_[facei];
1413 
1414  forAll(faces, i)
1415  {
1416  cop(result[facei], facei, fld[faces[i]], weights[i]);
1417  }
1418  }
1419  }
1420  }
1421 }
1422 
1423 
1424 template<class SourcePatch, class TargetPatch>
1425 template<class Type, class CombineOp>
1429  const Field<Type>& fld,
1430  const CombineOp& cop,
1431  const UList<Type>& defaultValues
1432 ) const
1433 {
1434  tmp<Field<Type>> tresult
1435  (
1436  new Field<Type>
1437  (
1438  srcAddress_.size(),
1439  Zero
1440  )
1441  );
1442 
1443  interpolateToSource
1444  (
1445  fld,
1447  tresult.ref(),
1448  defaultValues
1449  );
1450 
1451  return tresult;
1452 }
1453 
1454 
1455 template<class SourcePatch, class TargetPatch>
1456 template<class Type, class CombineOp>
1460  const tmp<Field<Type>>& tFld,
1461  const CombineOp& cop,
1462  const UList<Type>& defaultValues
1463 ) const
1464 {
1465  return interpolateToSource(tFld(), cop, defaultValues);
1466 }
1467 
1468 
1469 template<class SourcePatch, class TargetPatch>
1470 template<class Type, class CombineOp>
1474  const Field<Type>& fld,
1475  const CombineOp& cop,
1476  const UList<Type>& defaultValues
1477 ) const
1478 {
1479  tmp<Field<Type>> tresult
1480  (
1481  new Field<Type>
1482  (
1483  tgtAddress_.size(),
1484  Zero
1485  )
1486  );
1487 
1488  interpolateToTarget
1489  (
1490  fld,
1492  tresult.ref(),
1493  defaultValues
1494  );
1495 
1496  return tresult;
1497 }
1498 
1499 
1500 template<class SourcePatch, class TargetPatch>
1501 template<class Type, class CombineOp>
1505  const tmp<Field<Type>>& tFld,
1506  const CombineOp& cop,
1507  const UList<Type>& defaultValues
1508 ) const
1509 {
1510  return interpolateToTarget(tFld(), cop, defaultValues);
1511 }
1512 
1513 
1514 template<class SourcePatch, class TargetPatch>
1515 template<class Type>
1519  const Field<Type>& fld,
1520  const UList<Type>& defaultValues
1521 ) const
1522 {
1523  return interpolateToSource(fld, plusEqOp<Type>(), defaultValues);
1524 }
1525 
1526 
1527 template<class SourcePatch, class TargetPatch>
1528 template<class Type>
1532  const tmp<Field<Type>>& tFld,
1533  const UList<Type>& defaultValues
1534 ) const
1535 {
1536  return interpolateToSource(tFld(), plusEqOp<Type>(), defaultValues);
1537 }
1538 
1539 
1540 template<class SourcePatch, class TargetPatch>
1541 template<class Type>
1545  const Field<Type>& fld,
1546  const UList<Type>& defaultValues
1547 ) const
1548 {
1549  return interpolateToTarget(fld, plusEqOp<Type>(), defaultValues);
1550 }
1551 
1552 
1553 template<class SourcePatch, class TargetPatch>
1554 template<class Type>
1558  const tmp<Field<Type>>& tFld,
1559  const UList<Type>& defaultValues
1560 ) const
1561 {
1562  return interpolateToTarget(tFld(), plusEqOp<Type>(), defaultValues);
1563 }
1564 
1565 
1566 template<class SourcePatch, class TargetPatch>
1569  const SourcePatch& srcPatch,
1570  const TargetPatch& tgtPatch,
1571  const vector& n,
1572  const label tgtFacei,
1573  point& tgtPoint
1574 )
1575 const
1576 {
1577  const pointField& srcPoints = srcPatch.points();
1578 
1579  // Source face addresses that intersect target face tgtFacei
1580  const labelList& addr = tgtAddress_[tgtFacei];
1581 
1582  pointHit nearest;
1583  nearest.setDistance(great);
1584  label nearestFacei = -1;
1585 
1586  forAll(addr, i)
1587  {
1588  const label srcFacei = addr[i];
1589  const face& f = srcPatch[srcFacei];
1590 
1591  pointHit ray = f.ray(tgtPoint, n, srcPoints);
1592 
1593  if (ray.hit())
1594  {
1595  // tgtPoint = ray.rawPoint();
1596  return srcFacei;
1597  }
1598  else if (ray.distance() < nearest.distance())
1599  {
1600  nearest = ray;
1601  nearestFacei = srcFacei;
1602  }
1603  }
1604 
1605  if (nearest.hit() || nearest.eligibleMiss())
1606  {
1607  // tgtPoint = nearest.rawPoint();
1608  return nearestFacei;
1609  }
1610 
1611  return -1;
1612 }
1613 
1614 
1615 template<class SourcePatch, class TargetPatch>
1618  const SourcePatch& srcPatch,
1619  const TargetPatch& tgtPatch,
1620  const vector& n,
1621  const label srcFacei,
1622  point& srcPoint
1623 )
1624 const
1625 {
1626  const pointField& tgtPoints = tgtPatch.points();
1627 
1628  pointHit nearest;
1629  nearest.setDistance(great);
1630  label nearestFacei = -1;
1631 
1632  // Target face addresses that intersect source face srcFacei
1633  const labelList& addr = srcAddress_[srcFacei];
1634 
1635  forAll(addr, i)
1636  {
1637  const label tgtFacei = addr[i];
1638  const face& f = tgtPatch[tgtFacei];
1639 
1640  pointHit ray = f.ray(srcPoint, n, tgtPoints);
1641 
1642  if (ray.hit())
1643  {
1644  // srcPoint = ray.rawPoint();
1645  return tgtFacei;
1646  }
1647  else if (ray.distance() < nearest.distance())
1648  {
1649  nearest = ray;
1650  nearestFacei = tgtFacei;
1651  }
1652  }
1653 
1654  if (nearest.hit() || nearest.eligibleMiss())
1655  {
1656  // srcPoint = nearest.rawPoint();
1657  return nearestFacei;
1658  }
1659 
1660  return -1;
1661 }
1662 
1663 
1664 template<class SourcePatch, class TargetPatch>
1667  const SourcePatch& srcPatch,
1668  const TargetPatch& tgtPatch,
1669  const labelListList& srcAddress
1670 )
1671 const
1672 {
1673  OFstream os("faceConnectivity" + Foam::name(Pstream::myProcNo()) + ".obj");
1674 
1675  label ptI = 1;
1676 
1677  forAll(srcAddress, i)
1678  {
1679  const labelList& addr = srcAddress[i];
1680  const point& srcPt = srcPatch.faceCentres()[i];
1681 
1682  forAll(addr, j)
1683  {
1684  label tgtPtI = addr[j];
1685  const point& tgtPt = tgtPatch.faceCentres()[tgtPtI];
1686 
1687  meshTools::writeOBJ(os, srcPt);
1688  meshTools::writeOBJ(os, tgtPt);
1689 
1690  os << "l " << ptI << " " << ptI + 1 << endl;
1691 
1692  ptI += 2;
1693  }
1694  }
1695 }
1696 
1697 
1698 // ************************************************************************* //
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:428
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:226
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
Type gMin(const FieldField< Field, Type > &f)
Xfer< List< T > > xfer()
Transfer contents to the Xfer container.
Definition: ListI.H:177
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Output to file stream.
Definition: OFstream.H:82
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:163
static interpolationMethod wordTointerpolationMethod(const word &method)
Convert word to interpolationMethod.
const labelListList & tgtAddress() const
Return const access to target patch addressing.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
label tgtPointFace(const SourcePatch &srcPatch, const TargetPatch &tgtPatch, const vector &n, const label srcFacei, point &srcPoint) const
Return target patch face index of point on source patch face.
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.
label srcPointFace(const SourcePatch &srcPatch, const TargetPatch &tgtPatch, const vector &n, const label tgtFacei, point &tgtPoint) const
Return source patch face index of point on target patch face.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
interpolationMethod
Enumeration specifying interpolation method.
void writeFaceConnectivity(const SourcePatch &srcPatch, const TargetPatch &tgtPatch, const labelListList &srcAddress) const
Write face connectivity as OBJ file.
const scalarField & srcMagSf() const
Return const access to source patch face areas.
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
pointHit ray(const point &p, const vector &n, const pointField &, const intersection::algorithm alg=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return potential intersection with face with a ray starting.
void setDistance(const scalar d)
Definition: PointHit.H:186
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
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:184
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
bool hit() const
Is there a hit.
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
~AMIInterpolation()
Destructor.
void interpolateToTarget(const UList< Type > &fld, const CombineOp &cop, List< Type > &result, const UList< Type > &defaultValues=UList< Type >::null()) const
Interpolate from source to target with supplied op.
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
Base class for Arbitrary Mesh Interface (AMI) methods.
Definition: AMIMethod.H:55
volScalarField scalarField(fieldObject, mesh)
const scalarListList & srcWeights() const
Return const access to source patch weights.
static const char nl
Definition: Ostream.H:265
Type gMax(const FieldField< Field, Type > &f)
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
void interpolateToSource(const UList< Type > &fld, const CombineOp &cop, List< Type > &result, const UList< Type > &defaultValues=UList< Type >::null()) const
Interpolate from target to source with supplied op.
label constructSize() const
Constructed data size.
void setSize(const label)
Reset size of List.
Definition: List.C:281
Class containing processor-to-processor mapping information.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
Input from memory buffer stream.
Definition: IStringStream.H:49
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
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
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
void update(const SourcePatch &srcPatch, const TargetPatch &tgtPatch, const bool report)
Update addressing and weights.
bool eligibleMiss() const
Is this an eligible miss.
Definition: PointHit.H:164
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
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.
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
static tmp< scalarField > patchMagSf(const Patch &patch, const faceAreaIntersect::triangulationMode triMode)
Calculate the patch face magnitudes for the given tri-mode.