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