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-2018 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 "AMIMethod.H"
31 #include "meshTools.H"
32 #include "mapDistribute.H"
33 #include "flipOp.H"
34 #include "profiling.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 template<class SourcePatch, class TargetPatch>
39 const Foam::Enum
40 <
42  interpolationMethod
43 >
45 ({
46  { interpolationMethod::imDirect, "directAMI" },
47  { interpolationMethod::imMapNearest, "mapNearestAMI" },
48  { interpolationMethod::imFaceAreaWeight, "faceAreaWeightAMI" },
49  { interpolationMethod::imPartialFaceAreaWeight, "partialFaceAreaWeightAMI" }
50 });
51 
52 template<class SourcePatch, class TargetPatch>
54  false;
55 
56 template<class SourcePatch, class TargetPatch>
57 template<class Patch>
60 (
61  const Patch& patch,
63 )
64 {
65  tmp<scalarField> tResult(new scalarField(patch.size(), Zero));
66  scalarField& result = tResult.ref();
67 
68  const pointField& patchPoints = patch.localPoints();
69 
70  faceList patchFaceTris;
71 
72  forAll(result, patchFacei)
73  {
74  faceAreaIntersect::triangulate
75  (
76  patch.localFaces()[patchFacei],
77  patchPoints,
78  triMode,
79  patchFaceTris
80  );
81 
82  forAll(patchFaceTris, i)
83  {
84  result[patchFacei] +=
86  (
87  patchPoints[patchFaceTris[i][0]],
88  patchPoints[patchFaceTris[i][1]],
89  patchPoints[patchFaceTris[i][2]]
90  ).mag();
91  }
92  }
93 
94  return tResult;
95 }
96 
97 
98 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
99 
100 template<class SourcePatch, class TargetPatch>
102 (
103  const searchableSurface& surf,
104  pointField& pts
105 ) const
106 {
107  addProfiling(ami, "AMIInterpolation::projectPointsToSurface");
108 
109  DebugInfo<< "AMI: projecting points to surface" << endl;
110 
112 
113  surf.findNearest(pts, scalarField(pts.size(), GREAT), nearInfo);
114 
115  label nMiss = 0;
116  forAll(nearInfo, i)
117  {
118  const pointIndexHit& pi = nearInfo[i];
119 
120  if (pi.hit())
121  {
122  pts[i] = pi.hitPoint();
123  }
124  else
125  {
126  // Point remains unchanged
127  ++nMiss;
128  }
129  }
130 
131  if (nMiss > 0)
132  {
134  << "Error projecting points to surface: "
135  << nMiss << " faces could not be determined"
136  << abort(FatalError);
137  }
138 }
139 
140 
141 template<class SourcePatch, class TargetPatch>
143 (
144  const scalarList& patchAreas,
145  const word& patchName,
146  const labelListList& addr,
147  scalarListList& wght,
148  scalarField& wghtSum,
149  const bool conformal,
150  const bool output,
151  const scalar lowWeightTol
152 )
153 {
154  addProfiling(ami, "AMIInterpolation::normaliseWeights");
155 
156  // Normalise the weights
157  wghtSum.setSize(wght.size(), 0.0);
158  label nLowWeight = 0;
159 
160  forAll(wght, facei)
161  {
162  scalarList& w = wght[facei];
163 
164  if (w.size())
165  {
166  scalar denom = patchAreas[facei];
167 
168  scalar s = sum(w);
169  scalar t = s/denom;
170 
171  if (conformal)
172  {
173  denom = s;
174  }
175 
176  forAll(w, i)
177  {
178  w[i] /= denom;
179  }
180 
181  wghtSum[facei] = t;
182 
183  if (t < lowWeightTol)
184  {
185  nLowWeight++;
186  }
187  }
188  else
189  {
190  wghtSum[facei] = 0;
191  }
192  }
193 
194 
195  if (output)
196  {
197  const label nFace = returnReduce(wght.size(), sumOp<label>());
198 
199  if (nFace)
200  {
201  Info<< indent
202  << "AMI: Patch " << patchName
203  << " sum(weights)"
204  << " min:" << gMin(wghtSum)
205  << " max:" << gMax(wghtSum)
206  << " average:" << gAverage(wghtSum) << nl;
207 
208  const label nLow = returnReduce(nLowWeight, sumOp<label>());
209 
210  if (nLow)
211  {
212  Info<< indent
213  << "AMI: Patch " << patchName
214  << " identified " << nLow
215  << " faces with weights less than " << lowWeightTol
216  << endl;
217  }
218  }
219  }
220 }
221 
222 
223 template<class SourcePatch, class TargetPatch>
225 (
226  const autoPtr<mapDistribute>& targetMapPtr,
227  const scalarList& fineSrcMagSf,
228  const labelListList& fineSrcAddress,
229  const scalarListList& fineSrcWeights,
230 
231  const labelList& sourceRestrictAddressing,
232  const labelList& targetRestrictAddressing,
233 
234  scalarList& srcMagSf,
235  labelListList& srcAddress,
236  scalarListList& srcWeights,
237  scalarField& srcWeightsSum,
238  autoPtr<mapDistribute>& tgtMap
239 )
240 {
241  addProfiling(ami, "AMIInterpolation::agglomerate");
242 
243  label sourceCoarseSize =
244  (
245  sourceRestrictAddressing.size()
246  ? max(sourceRestrictAddressing)+1
247  : 0
248  );
249 
250  label targetCoarseSize =
251  (
252  targetRestrictAddressing.size()
253  ? max(targetRestrictAddressing)+1
254  : 0
255  );
256 
257  // Agglomerate face areas
258  {
259  srcMagSf.setSize(sourceRestrictAddressing.size(), 0.0);
260  forAll(sourceRestrictAddressing, facei)
261  {
262  label coarseFacei = sourceRestrictAddressing[facei];
263  srcMagSf[coarseFacei] += fineSrcMagSf[facei];
264  }
265  }
266 
267 
268  // Agglomerate weights and indices
269  if (targetMapPtr.valid())
270  {
271  const mapDistribute& map = targetMapPtr();
272 
273  // Get all restriction addressing.
274  labelList allRestrict(targetRestrictAddressing);
275  map.distribute(allRestrict);
276 
277  // So now we have agglomeration of the target side in
278  // allRestrict:
279  // 0..size-1 : local agglomeration (= targetRestrictAddressing
280  // (but potentially permutated))
281  // size.. : agglomeration data from other processors
282 
283 
284  // The trickiness in this algorithm is finding out the compaction
285  // of the remote data (i.e. allocation of the coarse 'slots'). We could
286  // either send across the slot compaction maps or just make sure
287  // that we allocate the slots in exactly the same order on both sending
288  // and receiving side (e.g. if the submap is set up to send 4 items,
289  // the constructMap is also set up to receive 4 items.
290 
291 
292  // Short note about the various types of indices:
293  // - face indices : indices into the geometry.
294  // - coarse face indices : how the faces get agglomerated
295  // - transferred data : how mapDistribute sends/receives data,
296  // - slots : indices into data after distribution (e.g. stencil,
297  // srcAddress/tgtAddress). Note: for fully local addressing
298  // the slots are equal to face indices.
299  // A mapDistribute has:
300  // - a subMap : these are face indices
301  // - a constructMap : these are from 'transferred-date' to slots
302 
303  labelListList tgtSubMap(Pstream::nProcs());
304 
305  // Local subMap is just identity
306  {
307  tgtSubMap[Pstream::myProcNo()] = identity(targetCoarseSize);
308  }
309 
310  forAll(map.subMap(), proci)
311  {
312  if (proci != Pstream::myProcNo())
313  {
314  // Combine entries that point to the same coarse element.
315  // The important bit is to loop over the data (and hand out
316  // compact indices ) in 'transferred data' order. This
317  // guarantees that we're doing exactly the
318  // same on sending and receiving side - e.g. the fourth element
319  // in the subMap is the fourth element received in the
320  // constructMap
321 
322  const labelList& elems = map.subMap()[proci];
323  const labelList& elemsMap =
324  map.constructMap()[Pstream::myProcNo()];
325  labelList& newSubMap = tgtSubMap[proci];
326  newSubMap.setSize(elems.size());
327 
328  labelList oldToNew(targetCoarseSize, -1);
329  label newi = 0;
330 
331  for (const label elemi : elems)
332  {
333  label fineElem = elemsMap[elemi];
334  label coarseElem = allRestrict[fineElem];
335  if (oldToNew[coarseElem] == -1)
336  {
337  oldToNew[coarseElem] = newi;
338  newSubMap[newi] = coarseElem;
339  ++newi;
340  }
341  }
342  newSubMap.setSize(newi);
343  }
344  }
345 
346  // Reconstruct constructMap by combining entries. Note that order
347  // of handing out indices should be the same as loop above to compact
348  // the sending map
349 
350  labelListList tgtConstructMap(Pstream::nProcs());
351 
352  // Local constructMap is just identity
353  {
354  tgtConstructMap[Pstream::myProcNo()] = identity(targetCoarseSize);
355  }
356 
357  labelList tgtCompactMap(map.constructSize());
358 
359  {
360  // Note that in special cases (e.g. 'appending' two AMIs) the
361  // local size after distributing can be longer than the number
362  // of faces. I.e. it duplicates elements.
363  // Since we don't know this size instead we loop over all
364  // reachable elements (using the local constructMap)
365 
366  const labelList& elemsMap = map.constructMap()[Pstream::myProcNo()];
367  for (const label fineElem : elemsMap)
368  {
369  label coarseElem = allRestrict[fineElem];
370  tgtCompactMap[fineElem] = coarseElem;
371  }
372  }
373 
374  label compacti = targetCoarseSize;
375 
376  // Compact data from other processors
377  forAll(map.constructMap(), proci)
378  {
379  if (proci != Pstream::myProcNo())
380  {
381  // Combine entries that point to the same coarse element. All
382  // elements now are remote data so we cannot use any local
383  // data here - use allRestrict instead.
384  const labelList& elems = map.constructMap()[proci];
385 
386  labelList& newConstructMap = tgtConstructMap[proci];
387  newConstructMap.setSize(elems.size());
388 
389  if (elems.size())
390  {
391  // Get the maximum target coarse size for this set of
392  // received data.
393  label remoteTargetCoarseSize = labelMin;
394  for (const label elemi : elems)
395  {
396  remoteTargetCoarseSize = max
397  (
398  remoteTargetCoarseSize,
399  allRestrict[elemi]
400  );
401  }
402  remoteTargetCoarseSize += 1;
403 
404  // Combine locally data coming from proci
405  labelList oldToNew(remoteTargetCoarseSize, -1);
406  label newi = 0;
407 
408  for (const label fineElem : elems)
409  {
410  // fineElem now points to section from proci
411  label coarseElem = allRestrict[fineElem];
412  if (oldToNew[coarseElem] == -1)
413  {
414  oldToNew[coarseElem] = newi;
415  tgtCompactMap[fineElem] = compacti;
416  newConstructMap[newi] = compacti++;
417  ++newi;
418  }
419  else
420  {
421  // Get compact index
422  label compacti = oldToNew[coarseElem];
423  tgtCompactMap[fineElem] = newConstructMap[compacti];
424  }
425  }
426  newConstructMap.setSize(newi);
427  }
428  }
429  }
430 
431  srcAddress.setSize(sourceCoarseSize);
432  srcWeights.setSize(sourceCoarseSize);
433 
434  forAll(fineSrcAddress, facei)
435  {
436  // All the elements contributing to facei. Are slots in
437  // mapDistribute'd data.
438  const labelList& elems = fineSrcAddress[facei];
439  const scalarList& weights = fineSrcWeights[facei];
440  const scalar fineArea = fineSrcMagSf[facei];
441 
442  label coarseFacei = sourceRestrictAddressing[facei];
443 
444  labelList& newElems = srcAddress[coarseFacei];
445  scalarList& newWeights = srcWeights[coarseFacei];
446 
447  forAll(elems, i)
448  {
449  label elemi = elems[i];
450  label coarseElemi = tgtCompactMap[elemi];
451 
452  label index = newElems.find(coarseElemi);
453  if (index == -1)
454  {
455  newElems.append(coarseElemi);
456  newWeights.append(fineArea*weights[i]);
457  }
458  else
459  {
460  newWeights[index] += fineArea*weights[i];
461  }
462  }
463  }
464 
465  tgtMap.reset
466  (
467  new mapDistribute
468  (
469  compacti,
470  std::move(tgtSubMap),
471  std::move(tgtConstructMap)
472  )
473  );
474  }
475  else
476  {
477  srcAddress.setSize(sourceCoarseSize);
478  srcWeights.setSize(sourceCoarseSize);
479 
480  forAll(fineSrcAddress, facei)
481  {
482  // All the elements contributing to facei. Are slots in
483  // mapDistribute'd data.
484  const labelList& elems = fineSrcAddress[facei];
485  const scalarList& weights = fineSrcWeights[facei];
486  const scalar fineArea = fineSrcMagSf[facei];
487 
488  label coarseFacei = sourceRestrictAddressing[facei];
489 
490  labelList& newElems = srcAddress[coarseFacei];
491  scalarList& newWeights = srcWeights[coarseFacei];
492 
493  forAll(elems, i)
494  {
495  label elemi = elems[i];
496  label coarseElemi = targetRestrictAddressing[elemi];
497 
498  label index = newElems.find(coarseElemi);
499  if (index == -1)
500  {
501  newElems.append(coarseElemi);
502  newWeights.append(fineArea*weights[i]);
503  }
504  else
505  {
506  newWeights[index] += fineArea*weights[i];
507  }
508  }
509  }
510  }
511 
512  // Weights normalisation
513  normaliseWeights
514  (
515  srcMagSf,
516  "source",
517  srcAddress,
518  srcWeights,
519  srcWeightsSum,
520  true,
521  false,
522  -1
523  );
524 }
525 
526 
527 template<class SourcePatch, class TargetPatch>
529 (
530  const SourcePatch& srcPatch,
531  const TargetPatch& tgtPatch,
532  const autoPtr<searchableSurface>& surfPtr
533 )
534 {
535  if (surfPtr.valid())
536  {
537  // Create new patches for source and target
538  pointField srcPoints = srcPatch.points();
539  SourcePatch srcPatch0
540  (
541  SubList<face>
542  (
543  srcPatch,
544  srcPatch.size(),
545  0
546  ),
547  srcPoints
548  );
549 
550  if (debug)
551  {
552  OFstream os("amiSrcPoints.obj");
553  for (const point& pt : srcPoints)
554  {
555  meshTools::writeOBJ(os, pt);
556  }
557  }
558 
559  pointField tgtPoints = tgtPatch.points();
560  TargetPatch tgtPatch0
561  (
562  SubList<face>
563  (
564  tgtPatch,
565  tgtPatch.size(),
566  0
567  ),
568  tgtPoints
569  );
570 
571  if (debug)
572  {
573  OFstream os("amiTgtPoints.obj");
574  for (const point& pt : tgtPoints)
575  {
576  meshTools::writeOBJ(os, pt);
577  }
578  }
579 
580 
581  // Map source and target patches onto projection surface
582  projectPointsToSurface(surfPtr(), srcPoints);
583  projectPointsToSurface(surfPtr(), tgtPoints);
584 
585 
586  // Calculate AMI interpolation
587  update(srcPatch0, tgtPatch0);
588  }
589  else
590  {
591  update(srcPatch, tgtPatch);
592  }
593 }
594 
595 
596 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
597 
598 template<class SourcePatch, class TargetPatch>
600 (
601  const SourcePatch& srcPatch,
602  const TargetPatch& tgtPatch,
604  const bool requireMatch,
605  const interpolationMethod& method,
606  const scalar lowWeightCorrection,
607  const bool reverseTarget
608 )
609 :
610  methodName_(interpolationMethodNames_[method]),
611  reverseTarget_(reverseTarget),
612  requireMatch_(requireMatch),
613  singlePatchProc_(-999),
614  lowWeightCorrection_(lowWeightCorrection),
615  srcMagSf_(),
616  srcAddress_(),
617  srcWeights_(),
618  srcWeightsSum_(),
619  tgtMagSf_(),
620  tgtAddress_(),
621  tgtWeights_(),
622  tgtWeightsSum_(),
623  triMode_(triMode),
624  srcMapPtr_(nullptr),
625  tgtMapPtr_(nullptr)
626 {
627  update(srcPatch, tgtPatch);
628 }
629 
630 
631 template<class SourcePatch, class TargetPatch>
633 (
634  const SourcePatch& srcPatch,
635  const TargetPatch& tgtPatch,
637  const bool requireMatch,
638  const word& methodName,
639  const scalar lowWeightCorrection,
640  const bool reverseTarget
641 )
642 :
643  methodName_(methodName),
644  reverseTarget_(reverseTarget),
645  requireMatch_(requireMatch),
646  singlePatchProc_(-999),
647  lowWeightCorrection_(lowWeightCorrection),
648  srcMagSf_(),
649  srcAddress_(),
650  srcWeights_(),
651  srcWeightsSum_(),
652  tgtMagSf_(),
653  tgtAddress_(),
654  tgtWeights_(),
655  tgtWeightsSum_(),
656  triMode_(triMode),
657  srcMapPtr_(nullptr),
658  tgtMapPtr_(nullptr)
659 {
660  update(srcPatch, tgtPatch);
661 }
662 
663 
664 template<class SourcePatch, class TargetPatch>
666 (
667  const SourcePatch& srcPatch,
668  const TargetPatch& tgtPatch,
669  const autoPtr<searchableSurface>& surfPtr,
671  const bool requireMatch,
672  const interpolationMethod& method,
673  const scalar lowWeightCorrection,
674  const bool reverseTarget
675 )
676 :
677  methodName_(interpolationMethodNames_[method]),
678  reverseTarget_(reverseTarget),
679  requireMatch_(requireMatch),
680  singlePatchProc_(-999),
681  lowWeightCorrection_(lowWeightCorrection),
682  srcMagSf_(),
683  srcAddress_(),
684  srcWeights_(),
685  srcWeightsSum_(),
686  tgtMagSf_(),
687  tgtAddress_(),
688  tgtWeights_(),
689  tgtWeightsSum_(),
690  triMode_(triMode),
691  srcMapPtr_(nullptr),
692  tgtMapPtr_(nullptr)
693 {
694  constructFromSurface(srcPatch, tgtPatch, surfPtr);
695 }
696 
697 
698 template<class SourcePatch, class TargetPatch>
700 (
701  const SourcePatch& srcPatch,
702  const TargetPatch& tgtPatch,
703  const autoPtr<searchableSurface>& surfPtr,
705  const bool requireMatch,
706  const word& methodName,
707  const scalar lowWeightCorrection,
708  const bool reverseTarget
709 )
710 :
711  methodName_(methodName),
712  reverseTarget_(reverseTarget),
713  requireMatch_(requireMatch),
714  singlePatchProc_(-999),
715  lowWeightCorrection_(lowWeightCorrection),
716  srcMagSf_(),
717  srcAddress_(),
718  srcWeights_(),
719  srcWeightsSum_(),
720  tgtMagSf_(),
721  tgtAddress_(),
722  tgtWeights_(),
723  tgtWeightsSum_(),
724  triMode_(triMode),
725  srcMapPtr_(nullptr),
726  tgtMapPtr_(nullptr)
727 {
728  constructFromSurface(srcPatch, tgtPatch, surfPtr);
729 }
730 
731 
732 template<class SourcePatch, class TargetPatch>
734 (
736  const labelList& sourceRestrictAddressing,
737  const labelList& targetRestrictAddressing
738 )
739 :
740  methodName_(fineAMI.methodName_),
741  reverseTarget_(fineAMI.reverseTarget_),
742  requireMatch_(fineAMI.requireMatch_),
743  singlePatchProc_(fineAMI.singlePatchProc_),
744  lowWeightCorrection_(-1.0),
745  srcMagSf_(),
746  srcAddress_(),
747  srcWeights_(),
748  srcWeightsSum_(),
749  tgtMagSf_(),
750  tgtAddress_(),
751  tgtWeights_(),
752  tgtWeightsSum_(),
753  triMode_(fineAMI.triMode_),
754  srcMapPtr_(nullptr),
755  tgtMapPtr_(nullptr)
756 {
757  label sourceCoarseSize =
758  (
759  sourceRestrictAddressing.size()
760  ? max(sourceRestrictAddressing)+1
761  : 0
762  );
763 
764  label neighbourCoarseSize =
765  (
766  targetRestrictAddressing.size()
767  ? max(targetRestrictAddressing)+1
768  : 0
769  );
770 
771  if (debug & 2)
772  {
773  Pout<< "AMI: Creating addressing and weights as agglomeration of AMI :"
774  << " source:" << fineAMI.srcAddress().size()
775  << " target:" << fineAMI.tgtAddress().size()
776  << " coarse source size:" << sourceCoarseSize
777  << " neighbour source size:" << neighbourCoarseSize
778  << endl;
779  }
780 
781  if
782  (
783  fineAMI.srcAddress().size() != sourceRestrictAddressing.size()
784  || fineAMI.tgtAddress().size() != targetRestrictAddressing.size()
785  )
786  {
788  << "Size mismatch." << nl
789  << "Source patch size:" << fineAMI.srcAddress().size() << nl
790  << "Source agglomeration size:"
791  << sourceRestrictAddressing.size() << nl
792  << "Target patch size:" << fineAMI.tgtAddress().size() << nl
793  << "Target agglomeration size:"
794  << targetRestrictAddressing.size()
795  << exit(FatalError);
796  }
797 
798 
799  // Agglomerate addresses and weights
800 
801  agglomerate
802  (
803  fineAMI.tgtMapPtr_,
804  fineAMI.srcMagSf(),
805  fineAMI.srcAddress(),
806  fineAMI.srcWeights(),
807 
808  sourceRestrictAddressing,
809  targetRestrictAddressing,
810 
811  srcMagSf_,
812  srcAddress_,
813  srcWeights_,
814  srcWeightsSum_,
815  tgtMapPtr_
816  );
817 
818  agglomerate
819  (
820  fineAMI.srcMapPtr_,
821  fineAMI.tgtMagSf(),
822  fineAMI.tgtAddress(),
823  fineAMI.tgtWeights(),
824 
825  targetRestrictAddressing,
826  sourceRestrictAddressing,
827 
828  tgtMagSf_,
829  tgtAddress_,
830  tgtWeights_,
831  tgtWeightsSum_,
832  srcMapPtr_
833  );
834 }
835 
836 
837 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
838 
839 template<class SourcePatch, class TargetPatch>
841 {}
842 
843 
844 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
845 
846 template<class SourcePatch, class TargetPatch>
848 (
849  const SourcePatch& srcPatch,
850  const TargetPatch& tgtPatch
851 )
852 {
853  addProfiling(ami, "AMIInterpolation::update");
854 
855  label srcTotalSize = returnReduce(srcPatch.size(), sumOp<label>());
856  label tgtTotalSize = returnReduce(tgtPatch.size(), sumOp<label>());
857 
858  if (srcTotalSize == 0)
859  {
860  DebugInfo<< "AMI: no source faces present - no addressing constructed"
861  << endl;
862 
863  return;
864  }
865 
866  Info<< indent
867  << "AMI: Creating addressing and weights between "
868  << srcTotalSize << " source faces and "
869  << tgtTotalSize << " target faces"
870  << endl;
871 
872  // Calculate if patches present on multiple processors
873  singlePatchProc_ = calcDistribution(srcPatch, tgtPatch);
874 
875  if (singlePatchProc_ == -1)
876  {
877  // Convert local addressing to global addressing
878  globalIndex globalSrcFaces(srcPatch.size());
879  globalIndex globalTgtFaces(tgtPatch.size());
880 
881  // Create processor map of overlapping faces. This map gets
882  // (possibly remote) faces from the tgtPatch such that they (together)
883  // cover all of the srcPatch
884  autoPtr<mapDistribute> mapPtr = calcProcMap(srcPatch, tgtPatch);
885  const mapDistribute& map = mapPtr();
886 
887  // Create new target patch that fully encompasses source patch
888 
889  // Faces and points
890  faceList newTgtFaces;
891  pointField newTgtPoints;
892 
893  // Original faces from tgtPatch (in globalIndexing since might be
894  // remote)
895  labelList tgtFaceIDs;
896  distributeAndMergePatches
897  (
898  map,
899  tgtPatch,
900  globalTgtFaces,
901  newTgtFaces,
902  newTgtPoints,
903  tgtFaceIDs
904  );
905 
906  const TargetPatch
907  newTgtPatch
908  (
910  (
911  newTgtFaces,
912  newTgtFaces.size()
913  ),
914  newTgtPoints
915  );
916 
917  // Calculate AMI interpolation
919  (
921  (
922  methodName_,
923  srcPatch,
924  newTgtPatch,
925  triMode_,
926  reverseTarget_,
927  requireMatch_ && (lowWeightCorrection_ < 0)
928  )
929  );
930 
931  AMIPtr->calculate
932  (
933  srcAddress_,
934  srcWeights_,
935  tgtAddress_,
936  tgtWeights_
937  );
938 
939 
940  // Note: using patch face areas calculated by the AMI method
941  // - TODO: move into the calculate or normalise method?
942  AMIPtr->setMagSf(tgtPatch, map, srcMagSf_, tgtMagSf_);
943 
944 
945  // Now
946  // ~~~
947  // srcAddress_ : per srcPatch face a list of the newTgtPatch (not
948  // tgtPatch) faces it overlaps
949  // tgtAddress_ : per newTgtPatch (not tgtPatch) face a list of the
950  // srcPatch faces it overlaps
951 
952  if (debug)
953  {
954  writeFaceConnectivity(srcPatch, newTgtPatch, srcAddress_);
955  }
956 
957 
958  // Rework newTgtPatch indices into globalIndices of tgtPatch
959  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
960 
961 
962  for (labelList& addressing : srcAddress_)
963  {
964  for (label& addr : addressing)
965  {
966  addr = tgtFaceIDs[addr];
967  }
968  }
969 
970  for (labelList& addressing : tgtAddress_)
971  {
972  globalSrcFaces.inplaceToGlobal(addressing);
973  }
974 
975  // Send data back to originating procs. Note that contributions
976  // from different processors get added (ListOps::appendEqOp)
977 
978  mapDistributeBase::distribute
979  (
980  Pstream::commsTypes::nonBlocking,
981  List<labelPair>(),
982  tgtPatch.size(),
983  map.constructMap(),
984  false, // has flip
985  map.subMap(),
986  false, // has flip
987  tgtAddress_,
989  flipOp(), // flip operation
990  labelList()
991  );
992 
993  mapDistributeBase::distribute
994  (
995  Pstream::commsTypes::nonBlocking,
996  List<labelPair>(),
997  tgtPatch.size(),
998  map.constructMap(),
999  false,
1000  map.subMap(),
1001  false,
1002  tgtWeights_,
1004  flipOp(),
1005  scalarList()
1006  );
1007 
1008  // weights normalisation
1009  AMIPtr->normaliseWeights(true, *this);
1010 
1011  // Cache maps and reset addresses
1012  List<Map<label>> cMap;
1013  srcMapPtr_.reset(new mapDistribute(globalSrcFaces, tgtAddress_, cMap));
1014  tgtMapPtr_.reset(new mapDistribute(globalTgtFaces, srcAddress_, cMap));
1015  }
1016  else
1017  {
1018  // Calculate AMI interpolation
1020  (
1022  (
1023  methodName_,
1024  srcPatch,
1025  tgtPatch,
1026  triMode_,
1027  reverseTarget_,
1028  requireMatch_ && (lowWeightCorrection_ < 0)
1029  )
1030  );
1031 
1032  AMIPtr->calculate
1033  (
1034  srcAddress_,
1035  srcWeights_,
1036  tgtAddress_,
1037  tgtWeights_
1038  );
1039 
1040  srcMagSf_.transfer(AMIPtr->srcMagSf());
1041  tgtMagSf_.transfer(AMIPtr->tgtMagSf());
1042 
1043  AMIPtr->normaliseWeights(true, *this);
1044  }
1045 
1046  if (debug)
1047  {
1048  Info<< "AMIInterpolation : Constructed addressing and weights" << nl
1049  << " triMode :"
1050  << faceAreaIntersect::triangulationModeNames_[triMode_] << nl
1051  << " singlePatchProc:" << singlePatchProc_ << nl
1052  << " srcMagSf :" << gSum(srcMagSf_) << nl
1053  << " tgtMagSf :" << gSum(tgtMagSf_) << nl
1054  << endl;
1055  }
1056 }
1057 
1058 
1059 template<class SourcePatch, class TargetPatch>
1062  const SourcePatch& srcPatch,
1063  const TargetPatch& tgtPatch
1064 )
1065 {
1066  addProfiling(ami, "AMIInterpolation::append");
1067 
1068  // Create a new interpolation
1070  (
1072  (
1073  srcPatch,
1074  tgtPatch,
1075  triMode_,
1076  requireMatch_,
1077  methodName_,
1078  lowWeightCorrection_,
1079  reverseTarget_
1080  )
1081  );
1082 
1083  // If parallel then combine the mapDistribution and re-index
1084  if (singlePatchProc_ == -1)
1085  {
1086  labelListList& srcSubMap = srcMapPtr_->subMap();
1087  labelListList& srcConstructMap = srcMapPtr_->constructMap();
1088 
1089  labelListList& tgtSubMap = tgtMapPtr_->subMap();
1090  labelListList& tgtConstructMap = tgtMapPtr_->constructMap();
1091 
1092  labelListList& newSrcSubMap = newPtr->srcMapPtr_->subMap();
1093  labelListList& newSrcConstructMap = newPtr->srcMapPtr_->constructMap();
1094 
1095  labelListList& newTgtSubMap = newPtr->tgtMapPtr_->subMap();
1096  labelListList& newTgtConstructMap = newPtr->tgtMapPtr_->constructMap();
1097 
1098  // Re-calculate the source indices
1099  {
1100  labelList mapMap(0), newMapMap(0);
1101  forAll(srcSubMap, proci)
1102  {
1103  mapMap.append
1104  (
1105  identity
1106  (
1107  srcConstructMap[proci].size(),
1108  (mapMap.size() + newMapMap.size())
1109  )
1110  );
1111  newMapMap.append
1112  (
1113  identity
1114  (
1115  newSrcConstructMap[proci].size(),
1116  (mapMap.size() + newMapMap.size())
1117  )
1118  );
1119  }
1120 
1121  forAll(srcSubMap, proci)
1122  {
1123  forAll(srcConstructMap[proci], srci)
1124  {
1125  srcConstructMap[proci][srci] =
1126  mapMap[srcConstructMap[proci][srci]];
1127  }
1128  }
1129 
1130  forAll(srcSubMap, proci)
1131  {
1132  forAll(newSrcConstructMap[proci], srci)
1133  {
1134  newSrcConstructMap[proci][srci] =
1135  newMapMap[newSrcConstructMap[proci][srci]];
1136  }
1137  }
1138 
1139  forAll(tgtAddress_, tgti)
1140  {
1141  forAll(tgtAddress_[tgti], tgtj)
1142  {
1143  tgtAddress_[tgti][tgtj] =
1144  mapMap[tgtAddress_[tgti][tgtj]];
1145  }
1146  }
1147 
1148  forAll(newPtr->tgtAddress_, tgti)
1149  {
1150  forAll(newPtr->tgtAddress_[tgti], tgtj)
1151  {
1152  newPtr->tgtAddress_[tgti][tgtj] =
1153  newMapMap[newPtr->tgtAddress_[tgti][tgtj]];
1154  }
1155  }
1156  }
1157 
1158  // Re-calculate the target indices
1159  {
1160  labelList mapMap(0), newMapMap(0);
1161  forAll(srcSubMap, proci)
1162  {
1163  mapMap.append
1164  (
1165  identity
1166  (
1167  tgtConstructMap[proci].size(),
1168  (mapMap.size() + newMapMap.size())
1169  )
1170  );
1171  newMapMap.append
1172  (
1173  identity
1174  (
1175  newTgtConstructMap[proci].size(),
1176  (mapMap.size() + newMapMap.size())
1177  )
1178  );
1179  }
1180 
1181  forAll(srcSubMap, proci)
1182  {
1183  forAll(tgtConstructMap[proci], tgti)
1184  {
1185  tgtConstructMap[proci][tgti] =
1186  mapMap[tgtConstructMap[proci][tgti]];
1187  }
1188  }
1189 
1190  forAll(srcSubMap, proci)
1191  {
1192  forAll(newTgtConstructMap[proci], tgti)
1193  {
1194  newTgtConstructMap[proci][tgti] =
1195  newMapMap[newTgtConstructMap[proci][tgti]];
1196  }
1197  }
1198 
1199  forAll(srcAddress_, srci)
1200  {
1201  forAll(srcAddress_[srci], srcj)
1202  {
1203  srcAddress_[srci][srcj] =
1204  mapMap[srcAddress_[srci][srcj]];
1205  }
1206  }
1207 
1208  forAll(newPtr->srcAddress_, srci)
1209  {
1210  forAll(newPtr->srcAddress_[srci], srcj)
1211  {
1212  newPtr->srcAddress_[srci][srcj] =
1213  newMapMap[newPtr->srcAddress_[srci][srcj]];
1214  }
1215  }
1216  }
1217 
1218  // Sum the construction sizes
1219  srcMapPtr_->constructSize() += newPtr->srcMapPtr_->constructSize();
1220  tgtMapPtr_->constructSize() += newPtr->tgtMapPtr_->constructSize();
1221 
1222  // Combine the maps
1223  forAll(srcSubMap, proci)
1224  {
1225  srcSubMap[proci].append(newSrcSubMap[proci]);
1226  srcConstructMap[proci].append(newSrcConstructMap[proci]);
1227 
1228  tgtSubMap[proci].append(newTgtSubMap[proci]);
1229  tgtConstructMap[proci].append(newTgtConstructMap[proci]);
1230  }
1231  }
1232 
1233  // Combine new and current source data
1234  forAll(srcMagSf_, srcFacei)
1235  {
1236  srcAddress_[srcFacei].append(newPtr->srcAddress()[srcFacei]);
1237  srcWeights_[srcFacei].append(newPtr->srcWeights()[srcFacei]);
1238  srcWeightsSum_[srcFacei] += newPtr->srcWeightsSum()[srcFacei];
1239  }
1240 
1241  // Combine new and current target data
1242  forAll(tgtMagSf_, tgtFacei)
1243  {
1244  tgtAddress_[tgtFacei].append(newPtr->tgtAddress()[tgtFacei]);
1245  tgtWeights_[tgtFacei].append(newPtr->tgtWeights()[tgtFacei]);
1246  tgtWeightsSum_[tgtFacei] += newPtr->tgtWeightsSum()[tgtFacei];
1247  }
1248 }
1249 
1250 
1251 template<class SourcePatch, class TargetPatch>
1254  const bool conformal,
1255  const bool output
1256 )
1257 {
1258  normaliseWeights
1259  (
1260  srcMagSf_,
1261  "source",
1262  srcAddress_,
1263  srcWeights_,
1264  srcWeightsSum_,
1265  conformal,
1266  output,
1267  lowWeightCorrection_
1268  );
1269 
1270  normaliseWeights
1271  (
1272  tgtMagSf_,
1273  "target",
1274  tgtAddress_,
1275  tgtWeights_,
1276  tgtWeightsSum_,
1277  conformal,
1278  output,
1279  lowWeightCorrection_
1280  );
1281 }
1282 
1283 
1284 template<class SourcePatch, class TargetPatch>
1285 template<class Type, class CombineOp>
1288  const UList<Type>& fld,
1289  const CombineOp& cop,
1290  List<Type>& result,
1291  const UList<Type>& defaultValues
1292 ) const
1293 {
1294  addProfiling(ami, "AMIInterpolation::interpolateToTarget");
1295 
1296  if (fld.size() != srcAddress_.size())
1297  {
1299  << "Supplied field size is not equal to source patch size" << nl
1300  << " source patch = " << srcAddress_.size() << nl
1301  << " target patch = " << tgtAddress_.size() << nl
1302  << " supplied field = " << fld.size()
1303  << abort(FatalError);
1304  }
1305 
1306  if (lowWeightCorrection_ > 0)
1307  {
1308  if (defaultValues.size() != tgtAddress_.size())
1309  {
1311  << "Employing default values when sum of weights falls below "
1312  << lowWeightCorrection_
1313  << " but supplied default field size is not equal to target "
1314  << "patch size" << nl
1315  << " default values = " << defaultValues.size() << nl
1316  << " target patch = " << tgtAddress_.size() << nl
1317  << abort(FatalError);
1318  }
1319  }
1320 
1321  result.setSize(tgtAddress_.size());
1322 
1323  if (singlePatchProc_ == -1)
1324  {
1325  const mapDistribute& map = srcMapPtr_();
1326 
1327  List<Type> work(fld);
1328  map.distribute(work);
1329 
1330  forAll(result, facei)
1331  {
1332  if (tgtWeightsSum_[facei] < lowWeightCorrection_)
1333  {
1334  result[facei] = defaultValues[facei];
1335  }
1336  else
1337  {
1338  const labelList& faces = tgtAddress_[facei];
1339  const scalarList& weights = tgtWeights_[facei];
1340 
1341  forAll(faces, i)
1342  {
1343  cop(result[facei], facei, work[faces[i]], weights[i]);
1344  }
1345  }
1346  }
1347  }
1348  else
1349  {
1350  forAll(result, facei)
1351  {
1352  if (tgtWeightsSum_[facei] < lowWeightCorrection_)
1353  {
1354  result[facei] = defaultValues[facei];
1355  }
1356  else
1357  {
1358  const labelList& faces = tgtAddress_[facei];
1359  const scalarList& weights = tgtWeights_[facei];
1360 
1361  forAll(faces, i)
1362  {
1363  cop(result[facei], facei, fld[faces[i]], weights[i]);
1364  }
1365  }
1366  }
1367  }
1368 }
1369 
1370 
1371 template<class SourcePatch, class TargetPatch>
1372 template<class Type, class CombineOp>
1375  const UList<Type>& fld,
1376  const CombineOp& cop,
1377  List<Type>& result,
1378  const UList<Type>& defaultValues
1379 ) const
1380 {
1381  addProfiling(ami, "AMIInterpolation::interpolateToSource");
1382 
1383  if (fld.size() != tgtAddress_.size())
1384  {
1386  << "Supplied field size is not equal to target patch size" << nl
1387  << " source patch = " << srcAddress_.size() << nl
1388  << " target patch = " << tgtAddress_.size() << nl
1389  << " supplied field = " << fld.size()
1390  << abort(FatalError);
1391  }
1392 
1393  if (lowWeightCorrection_ > 0)
1394  {
1395  if (defaultValues.size() != srcAddress_.size())
1396  {
1398  << "Employing default values when sum of weights falls below "
1399  << lowWeightCorrection_
1400  << " but supplied default field size is not equal to target "
1401  << "patch size" << nl
1402  << " default values = " << defaultValues.size() << nl
1403  << " source patch = " << srcAddress_.size() << nl
1404  << abort(FatalError);
1405  }
1406  }
1407 
1408  result.setSize(srcAddress_.size());
1409 
1410  if (singlePatchProc_ == -1)
1411  {
1412  const mapDistribute& map = tgtMapPtr_();
1413 
1414  List<Type> work(fld);
1415  map.distribute(work);
1416 
1417  forAll(result, facei)
1418  {
1419  if (srcWeightsSum_[facei] < lowWeightCorrection_)
1420  {
1421  result[facei] = defaultValues[facei];
1422  }
1423  else
1424  {
1425  const labelList& faces = srcAddress_[facei];
1426  const scalarList& weights = srcWeights_[facei];
1427 
1428  forAll(faces, i)
1429  {
1430  cop(result[facei], facei, work[faces[i]], weights[i]);
1431  }
1432  }
1433  }
1434  }
1435  else
1436  {
1437  forAll(result, facei)
1438  {
1439  if (srcWeightsSum_[facei] < lowWeightCorrection_)
1440  {
1441  result[facei] = defaultValues[facei];
1442  }
1443  else
1444  {
1445  const labelList& faces = srcAddress_[facei];
1446  const scalarList& weights = srcWeights_[facei];
1447 
1448  forAll(faces, i)
1449  {
1450  cop(result[facei], facei, fld[faces[i]], weights[i]);
1451  }
1452  }
1453  }
1454  }
1455 }
1456 
1457 
1458 template<class SourcePatch, class TargetPatch>
1459 template<class Type, class CombineOp>
1463  const Field<Type>& fld,
1464  const CombineOp& cop,
1465  const UList<Type>& defaultValues
1466 ) const
1467 {
1468  tmp<Field<Type>> tresult
1469  (
1470  new Field<Type>
1471  (
1472  srcAddress_.size(),
1473  Zero
1474  )
1475  );
1476 
1477  interpolateToSource
1478  (
1479  fld,
1481  tresult.ref(),
1482  defaultValues
1483  );
1484 
1485  return tresult;
1486 }
1487 
1488 
1489 template<class SourcePatch, class TargetPatch>
1490 template<class Type, class CombineOp>
1494  const tmp<Field<Type>>& tFld,
1495  const CombineOp& cop,
1496  const UList<Type>& defaultValues
1497 ) const
1498 {
1499  return interpolateToSource(tFld(), cop, defaultValues);
1500 }
1501 
1502 
1503 template<class SourcePatch, class TargetPatch>
1504 template<class Type, class CombineOp>
1508  const Field<Type>& fld,
1509  const CombineOp& cop,
1510  const UList<Type>& defaultValues
1511 ) const
1512 {
1513  tmp<Field<Type>> tresult
1514  (
1515  new Field<Type>
1516  (
1517  tgtAddress_.size(),
1518  Zero
1519  )
1520  );
1521 
1522  interpolateToTarget
1523  (
1524  fld,
1526  tresult.ref(),
1527  defaultValues
1528  );
1529 
1530  return tresult;
1531 }
1532 
1533 
1534 template<class SourcePatch, class TargetPatch>
1535 template<class Type, class CombineOp>
1539  const tmp<Field<Type>>& tFld,
1540  const CombineOp& cop,
1541  const UList<Type>& defaultValues
1542 ) const
1543 {
1544  return interpolateToTarget(tFld(), cop, defaultValues);
1545 }
1546 
1547 
1548 template<class SourcePatch, class TargetPatch>
1549 template<class Type>
1553  const Field<Type>& fld,
1554  const UList<Type>& defaultValues
1555 ) const
1556 {
1557  return interpolateToSource(fld, plusEqOp<Type>(), defaultValues);
1558 }
1559 
1560 
1561 template<class SourcePatch, class TargetPatch>
1562 template<class Type>
1566  const tmp<Field<Type>>& tFld,
1567  const UList<Type>& defaultValues
1568 ) const
1569 {
1570  return interpolateToSource(tFld(), plusEqOp<Type>(), defaultValues);
1571 }
1572 
1573 
1574 template<class SourcePatch, class TargetPatch>
1575 template<class Type>
1579  const Field<Type>& fld,
1580  const UList<Type>& defaultValues
1581 ) const
1582 {
1583  return interpolateToTarget(fld, plusEqOp<Type>(), defaultValues);
1584 }
1585 
1586 
1587 template<class SourcePatch, class TargetPatch>
1588 template<class Type>
1592  const tmp<Field<Type>>& tFld,
1593  const UList<Type>& defaultValues
1594 ) const
1595 {
1596  return interpolateToTarget(tFld(), plusEqOp<Type>(), defaultValues);
1597 }
1598 
1599 
1600 template<class SourcePatch, class TargetPatch>
1603  const SourcePatch& srcPatch,
1604  const TargetPatch& tgtPatch,
1605  const vector& n,
1606  const label tgtFacei,
1607  point& tgtPoint
1608 )
1609 const
1610 {
1611  const pointField& srcPoints = srcPatch.points();
1612 
1613  // Source face addresses that intersect target face tgtFacei
1614  const labelList& addr = tgtAddress_[tgtFacei];
1615 
1616  pointHit nearest;
1617  nearest.setDistance(GREAT);
1618  label nearestFacei = -1;
1619 
1620  for (const label srcFacei : addr)
1621  {
1622  const face& f = srcPatch[srcFacei];
1623 
1624  pointHit ray = f.ray(tgtPoint, n, srcPoints);
1625 
1626  if (ray.hit())
1627  {
1628  // tgtPoint = ray.rawPoint();
1629  return srcFacei;
1630  }
1631  else if (ray.distance() < nearest.distance())
1632  {
1633  nearest = ray;
1634  nearestFacei = srcFacei;
1635  }
1636  }
1637 
1638  if (nearest.hit() || nearest.eligibleMiss())
1639  {
1640  // tgtPoint = nearest.rawPoint();
1641  return nearestFacei;
1642  }
1643 
1644  return -1;
1645 }
1646 
1647 
1648 template<class SourcePatch, class TargetPatch>
1651  const SourcePatch& srcPatch,
1652  const TargetPatch& tgtPatch,
1653  const vector& n,
1654  const label srcFacei,
1655  point& srcPoint
1656 )
1657 const
1658 {
1659  const pointField& tgtPoints = tgtPatch.points();
1660 
1661  pointHit nearest;
1662  nearest.setDistance(GREAT);
1663  label nearestFacei = -1;
1664 
1665  // Target face addresses that intersect source face srcFacei
1666  const labelList& addr = srcAddress_[srcFacei];
1667 
1668  for (const label tgtFacei : addr)
1669  {
1670  const face& f = tgtPatch[tgtFacei];
1671 
1672  pointHit ray = f.ray(srcPoint, n, tgtPoints);
1673 
1674  if (ray.hit() || ray.eligibleMiss())
1675  {
1676  // srcPoint = ray.rawPoint();
1677  return tgtFacei;
1678  }
1679  else if (ray.distance() < nearest.distance())
1680  {
1681  nearest = ray;
1682  nearestFacei = tgtFacei;
1683  }
1684  }
1685 
1686  if (nearest.hit() || nearest.eligibleMiss())
1687  {
1688  // srcPoint = nearest.rawPoint();
1689  return nearestFacei;
1690  }
1691 
1692  return -1;
1693 }
1694 
1695 
1696 template<class SourcePatch, class TargetPatch>
1699  const SourcePatch& srcPatch,
1700  const TargetPatch& tgtPatch,
1701  const labelListList& srcAddress
1702 )
1703 const
1704 {
1705  OFstream os("faceConnectivity" + Foam::name(Pstream::myProcNo()) + ".obj");
1706 
1707  label pti = 1;
1708 
1709  forAll(srcAddress, i)
1710  {
1711  const labelList& addr = srcAddress[i];
1712  const point& srcPt = srcPatch.faceCentres()[i];
1713 
1714  for (const label tgtPti : addr)
1715  {
1716  const point& tgtPt = tgtPatch.faceCentres()[tgtPti];
1717 
1718  meshTools::writeOBJ(os, srcPt);
1719  meshTools::writeOBJ(os, tgtPt);
1720 
1721  os << "l " << pti << " " << pti + 1 << endl;
1722 
1723  pti += 2;
1724  }
1725  }
1726 }
1727 
1728 
1729 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::AMIInterpolation::interpolationMethod
interpolationMethod
Enumeration specifying interpolation method.
Definition: AMIInterpolation.H:90
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
profiling.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::scalarListList
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
Foam::mapDistributeBase::subMap
const labelListList & subMap() const
From subsetted data back to original data.
Definition: mapDistributeBase.H:282
Foam::AMIInterpolation::append
void append(const SourcePatch &srcPatch, const TargetPatch &tgtPatch)
Append additional addressing and weights.
Definition: AMIInterpolation.C:1061
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
meshTools.H
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::nearInfo
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
Definition: sampledTriSurfaceMesh.C:67
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::AMIInterpolation::interpolationMethodNames_
static const Enum< interpolationMethod > interpolationMethodNames_
Definition: AMIInterpolation.H:98
Foam::PointHit::hit
bool hit() const
Is there a hit.
Definition: PointHit.H:124
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
s
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))
Definition: gmvOutputSpray.H:25
Foam::AMIInterpolation::tgtPointFace
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.
Definition: AMIInterpolation.C:1650
update
mesh update()
Foam::faceAreaIntersect::triangulationMode
triangulationMode
Definition: faceAreaIntersect.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::PointHit< point >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
AMIInterpolation.H
Foam::AMIInterpolation::tgtAddress
const labelListList & tgtAddress() const
Return const access to target patch addressing.
Definition: AMIInterpolationI.H:144
Foam::AMIInterpolation::tgtWeights
const scalarListList & tgtWeights() const
Return const access to target patch weights.
Definition: AMIInterpolationI.H:160
Foam::AMIMethod
Base class for Arbitrary Mesh Interface (AMI) methods.
Definition: AMIMethod.H:60
Foam::AMIInterpolation::interpolateToTarget
void interpolateToTarget(const UList< Type > &fld, const CombineOp &cop, List< Type > &result, const UList< Type > &defaultValues=UList< Type >::null()) const
Definition: AMIInterpolation.C:1287
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:55
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:143
Foam::label
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:62
Foam::AMIInterpolation::interpolateToSource
void interpolateToSource(const UList< Type > &fld, const CombineOp &cop, List< Type > &result, const UList< Type > &defaultValues=UList< Type >::null()) const
Definition: AMIInterpolation.C:1374
Foam::Field< scalar >
Foam::PointHit::eligibleMiss
bool eligibleMiss() const
Is this an eligible miss.
Definition: PointHit.H:168
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::AMIInterpolation::~AMIInterpolation
~AMIInterpolation()
Destructor.
Definition: AMIInterpolation.C:840
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::AMIInterpolation::writeFaceConnectivity
void writeFaceConnectivity(const SourcePatch &srcPatch, const TargetPatch &tgtPatch, const labelListList &srcAddress) const
Write face connectivity as OBJ file.
Definition: AMIInterpolation.C:1698
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fld
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:152
addProfiling
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
Definition: profilingTrigger.H:114
Foam::FatalError
error FatalError
Foam::AMIInterpolation::patchMagSf
static tmp< scalarField > patchMagSf(const Patch &patch, const faceAreaIntersect::triangulationMode triMode)
Calculate the patch face magnitudes for the given tri-mode.
flipOp.H
Foam::searchableSurface::findNearest
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:307
Foam::AMIInterpolation::srcPointFace
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.
Definition: AMIInterpolation.C:1602
Foam::AMIInterpolation::srcMagSf
const List< scalar > & srcMagSf() const
Return const access to source patch face areas.
Definition: AMIInterpolationI.H:56
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:99
Foam::ListOps::appendEqOp
List helper to append y elements onto the end of x.
Definition: ListOps.H:569
Foam::multiplyWeightedOp
Definition: ops.H:250
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:372
mapDistribute.H
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::List< face >
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:81
Foam::UList< Type >
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
AMIMethod.H
Foam::labelMin
constexpr label labelMin
Definition: label.H:64
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::AMIInterpolation::update
void update(const SourcePatch &srcPatch, const TargetPatch &tgtPatch)
Update addressing and weights.
Definition: AMIInterpolation.C:848
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::plusEqOp
Definition: ops.H:72
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:57
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::AMIInterpolation::tgtMagSf
const List< scalar > & tgtMagSf() const
Return const access to target patch face areas.
Definition: AMIInterpolationI.H:128
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::mapDistributeBase::constructMap
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
Definition: mapDistributeBase.H:294
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::AMIInterpolation::srcAddress
const labelListList & srcAddress() const
Return const access to source patch addressing.
Definition: AMIInterpolationI.H:72
Foam::AMIInterpolation::srcWeights
const scalarListList & srcWeights() const
Return const access to source patch weights.
Definition: AMIInterpolationI.H:88
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:78
Foam::flipOp
Functor to negate primitives. Dummy for most other types.
Definition: flipOp.H:53
Foam::PointHit::setDistance
void setDistance(const scalar d)
Definition: PointHit.H:192