faceAreaWeightAMI.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2021 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 "faceAreaWeightAMI.H"
30 #include "profiling.H"
31 #include "OBJstream.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(faceAreaWeightAMI, 0);
39  addToRunTimeSelectionTable(AMIInterpolation, faceAreaWeightAMI, dict);
40  addToRunTimeSelectionTable(AMIInterpolation, faceAreaWeightAMI, component);
41 
42  // Backwards compatibility for pre v2106 versions
43  // - partialFaceAreaWeightAMI deprecated in v2106
45  (
46  AMIInterpolation,
47  faceAreaWeightAMI,
48  dict,
49  faceAreaWeightAMI,
50  partialFaceAreaWeightAMI,
51  2012
52  );
53 }
54 
55 
56 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
57 
58 /*
59  if (debug)
60  {
61  static label nAMI = 0;
62 
63  // Write out triangulated surfaces as OBJ files
64  OBJstream srcTriObj("srcTris_" + Foam::name(nAMI) + ".obj");
65  const pointField& srcPts = src.points();
66  forAll(srcTris_, facei)
67  {
68  const DynamicList<face>& faces = srcTris_[facei];
69  for (const face& f : faces)
70  {
71  srcTriObj.write
72  (
73  triPointRef(srcPts[f[0]], srcPts[f[1]], srcPts[f[2]])
74  );
75  }
76  }
77 
78  OBJstream tgtTriObj("tgtTris_" + Foam::name(nAMI) + ".obj");
79  const pointField& tgtPts = tgt.points();
80  forAll(tgtTris_, facei)
81  {
82  const DynamicList<face>& faces = tgtTris_[facei];
83  for (const face& f : faces)
84  {
85  tgtTriObj.write
86  (
87  triPointRef(tgtPts[f[0]], tgtPts[f[1]], tgtPts[f[2]])
88  );
89  }
90  }
91 
92  ++nAMI;
93  }
94 */
95 
96 
98 (
99  List<DynamicList<label>>& srcAddr,
100  List<DynamicList<scalar>>& srcWght,
101  List<DynamicList<point>>& srcCtr,
102  List<DynamicList<label>>& tgtAddr,
103  List<DynamicList<scalar>>& tgtWght,
104  label srcFacei,
105  label tgtFacei
106 )
107 {
108  addProfiling(ami, "faceAreaWeightAMI::calcAddressing");
109 
110  // Construct weights and addressing
111  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
112 
113  label nFacesRemaining = srcAddr.size();
114 
115  // List of tgt face neighbour faces
116  DynamicList<label> nbrFaces(10);
117 
118  // List of faces currently visited for srcFacei to avoid multiple hits
119  DynamicList<label> visitedFaces(10);
120 
121  // List to keep track of tgt faces used to seed src faces
122  labelList seedFaces(nFacesRemaining, -1);
123  seedFaces[srcFacei] = tgtFacei;
124 
125  // List to keep track of whether src face can be mapped
126  bitSet mapFlag(nFacesRemaining, true);
127 
128  // Reset starting seed
129  label startSeedi = 0;
130 
131  // Should all faces be matched?
132  const bool mustMatch = mustMatchFaces();
133 
134  bool continueWalk = true;
135  DynamicList<label> nonOverlapFaces;
136  do
137  {
138  nbrFaces.clear();
139  visitedFaces.clear();
140 
141  // Do advancing front starting from srcFacei,tgtFacei
142  bool faceProcessed = processSourceFace
143  (
144  srcFacei,
145  tgtFacei,
146 
147  nbrFaces,
148  visitedFaces,
149 
150  srcAddr,
151  srcWght,
152  srcCtr,
153  tgtAddr,
154  tgtWght
155  );
156 
157  mapFlag.unset(srcFacei);
158 
159  if (!faceProcessed)
160  {
161  nonOverlapFaces.append(srcFacei);
162  }
163 
164  // Choose new src face from current src face neighbour
165  continueWalk = setNextFaces
166  (
167  startSeedi,
168  srcFacei,
169  tgtFacei,
170  mapFlag,
171  seedFaces,
172  visitedFaces,
173  mustMatch
174  // pass in nonOverlapFaces for failed tree search?
175  );
176  } while (continueWalk);
177 
178  srcNonOverlap_.transfer(nonOverlapFaces);
179 }
180 
181 
183 (
184  const label srcFacei,
185  const label tgtStartFacei,
186 
187  // list of tgt face neighbour faces
188  DynamicList<label>& nbrFaces,
189  // list of faces currently visited for srcFacei to avoid multiple hits
190  DynamicList<label>& visitedFaces,
191 
192  // temporary storage for addressing, weights and centroid
193  List<DynamicList<label>>& srcAddr,
194  List<DynamicList<scalar>>& srcWght,
195  List<DynamicList<point>>& srcCtr,
196  List<DynamicList<label>>& tgtAddr,
197  List<DynamicList<scalar>>& tgtWght
198 )
199 {
200  addProfiling(ami, "faceAreaWeightAMI::processSourceFace");
201 
202  if (tgtStartFacei == -1)
203  {
204  return false;
205  }
206 
207  const auto& tgtPatch = this->tgtPatch();
208 
209  // append initial target face and neighbours
210  nbrFaces.append(tgtStartFacei);
211  appendNbrFaces(tgtStartFacei, tgtPatch, visitedFaces, nbrFaces);
212 
213  bool faceProcessed = false;
214 
215  label maxNeighbourFaces = nbrFaces.size();
216 
217  do
218  {
219  // process new target face
220  label tgtFacei = nbrFaces.remove();
221  visitedFaces.append(tgtFacei);
222 
223  scalar interArea = 0;
224  vector interCentroid(Zero);
225  calcInterArea(srcFacei, tgtFacei, interArea, interCentroid);
226 
227  // store when intersection fractional area > tolerance
228  if (interArea/srcMagSf_[srcFacei] > faceAreaIntersect::tolerance())
229  {
230  srcAddr[srcFacei].append(tgtFacei);
231  srcWght[srcFacei].append(interArea);
232  srcCtr[srcFacei].append(interCentroid);
233 
234  tgtAddr[tgtFacei].append(srcFacei);
235  tgtWght[tgtFacei].append(interArea);
236 
237  appendNbrFaces(tgtFacei, tgtPatch, visitedFaces, nbrFaces);
238 
239  faceProcessed = true;
240 
241  maxNeighbourFaces = max(maxNeighbourFaces, nbrFaces.size());
242  }
243 
244  } while (nbrFaces.size() > 0);
245 
246  if (debug > 1)
247  {
248  DebugVar(maxNeighbourFaces);
249  }
250 
251  return faceProcessed;
252 }
253 
254 
256 (
257  label& startSeedi,
258  label& srcFacei,
259  label& tgtFacei,
260  const bitSet& mapFlag,
261  labelList& seedFaces,
262  const DynamicList<label>& visitedFaces,
263  const bool errorOnNotFound
264 ) const
265 {
266  addProfiling(ami, "faceAreaWeightAMI::setNextFaces");
267 
268  if (mapFlag.count() == 0)
269  {
270  // No more faces to map
271  return false;
272  }
273 
274  const labelList& srcNbrFaces = this->srcPatch().faceFaces()[srcFacei];
275 
276  // Initialise tgtFacei
277  tgtFacei = -1;
278 
279  // Set possible seeds for later use
280  bool valuesSet = false;
281  for (label faceS: srcNbrFaces)
282  {
283  if (mapFlag.test(faceS) && seedFaces[faceS] == -1)
284  {
285  for (label faceT : visitedFaces)
286  {
287  const scalar threshold =
288  srcMagSf_[faceS]*faceAreaIntersect::tolerance();
289 
290  // Store when intersection fractional area > threshold
291  if (overlaps(faceS, faceT, threshold))
292  {
293  seedFaces[faceS] = faceT;
294 
295  if (!valuesSet)
296  {
297  srcFacei = faceS;
298  tgtFacei = faceT;
299  valuesSet = true;
300  }
301  }
302  }
303  }
304  }
305 
306  if (valuesSet)
307  {
308  return true;
309  }
310 
311  // Set next src and tgt faces if not set above
312  // - try to use existing seed
313  label facei = startSeedi;
314  if (!mapFlag.test(startSeedi))
315  {
316  facei = mapFlag.find_next(facei);
317  }
318  const label startSeedi0 = facei;
319 
320  bool foundNextSeed = false;
321  while (facei != -1)
322  {
323  if (!foundNextSeed)
324  {
325  startSeedi = facei;
326  foundNextSeed = true;
327  }
328 
329  if (seedFaces[facei] != -1)
330  {
331  srcFacei = facei;
332  tgtFacei = seedFaces[facei];
333 
334  return true;
335  }
336 
337  facei = mapFlag.find_next(facei);
338  }
339 
340  // Perform new search to find match
341  if (debug)
342  {
343  Pout<< "Advancing front stalled: searching for new "
344  << "target face" << endl;
345  }
346 
347  facei = startSeedi0;
348  while (facei != -1)
349  {
350  srcFacei = facei;
351  tgtFacei = findTargetFace(srcFacei, visitedFaces);
352 
353  if (tgtFacei >= 0)
354  {
355  return true;
356  }
357 
358  facei = mapFlag.find_next(facei);
359  }
360 
361  if (errorOnNotFound)
362  {
364  << "Unable to set target face for source face " << srcFacei
365  << abort(FatalError);
366  }
367 
368  return false;
369 }
370 
371 
373 (
374  const label srcFacei,
375  const label tgtFacei,
376  scalar& area,
377  vector& centroid
378 ) const
379 {
380  addProfiling(ami, "faceAreaWeightAMI::interArea");
381 
382  // Quick reject if either face has zero area
383  if
384  (
385  (srcMagSf_[srcFacei] < ROOTVSMALL)
386  || (tgtMagSf_[tgtFacei] < ROOTVSMALL)
387  )
388  {
389  return;
390  }
391 
392  const auto& srcPatch = this->srcPatch();
393  const auto& tgtPatch = this->tgtPatch();
394 
395  const pointField& srcPoints = srcPatch.points();
396  const pointField& tgtPoints = tgtPatch.points();
397 
398  // Create intersection object
399  faceAreaIntersect inter
400  (
401  srcPoints,
402  tgtPoints,
403  srcTris_[srcFacei],
404  tgtTris_[tgtFacei],
405  reverseTarget_,
406  AMIInterpolation::cacheIntersections_
407  );
408 
409  // Crude resultant norm
410  vector n(-srcPatch.faceNormals()[srcFacei]);
411  if (reverseTarget_)
412  {
413  n -= tgtPatch.faceNormals()[tgtFacei];
414  }
415  else
416  {
417  n += tgtPatch.faceNormals()[tgtFacei];
418  }
419  scalar magN = mag(n);
420 
421  const face& src = srcPatch[srcFacei];
422  const face& tgt = tgtPatch[tgtFacei];
423 
424  if (magN > ROOTVSMALL)
425  {
426  inter.calc(src, tgt, n/magN, area, centroid);
427  }
428  else
429  {
431  << "Invalid normal for source face " << srcFacei
432  << " points " << UIndirectList<point>(srcPoints, src)
433  << " target face " << tgtFacei
434  << " points " << UIndirectList<point>(tgtPoints, tgt)
435  << endl;
436  }
437 
438  if (AMIInterpolation::cacheIntersections_ && debug)
439  {
440  static OBJstream tris("intersectionTris.obj");
441  const auto& triPts = inter.triangles();
442  for (const auto& tp : triPts)
443  {
444  tris.write(triPointRef(tp[0], tp[1], tp[2]), false);
445  }
446  }
447 
448  if ((debug > 1) && (area > 0))
449  {
450  writeIntersectionOBJ(area, src, tgt, srcPoints, tgtPoints);
451  }
452 }
453 
454 
456 (
457  const label srcFacei,
458  const label tgtFacei,
459  const scalar threshold
460 ) const
461 {
462  // Quick reject if either face has zero area
463  if
464  (
465  (srcMagSf_[srcFacei] < ROOTVSMALL)
466  || (tgtMagSf_[tgtFacei] < ROOTVSMALL)
467  )
468  {
469  return false;
470  }
471 
472  const auto& srcPatch = this->srcPatch();
473  const auto& tgtPatch = this->tgtPatch();
474 
475  const pointField& srcPoints = srcPatch.points();
476  const pointField& tgtPoints = tgtPatch.points();
477 
478  faceAreaIntersect inter
479  (
480  srcPoints,
481  tgtPoints,
482  srcTris_[srcFacei],
483  tgtTris_[tgtFacei],
484  reverseTarget_,
485  AMIInterpolation::cacheIntersections_
486  );
487 
488  // Crude resultant norm
489  vector n(-srcPatch.faceNormals()[srcFacei]);
490  if (reverseTarget_)
491  {
492  n -= tgtPatch.faceNormals()[tgtFacei];
493  }
494  else
495  {
496  n += tgtPatch.faceNormals()[tgtFacei];
497  }
498  scalar magN = mag(n);
499 
500  const face& src = srcPatch[srcFacei];
501  const face& tgt = tgtPatch[tgtFacei];
502 
503  if (magN > ROOTVSMALL)
504  {
505  return inter.overlaps(src, tgt, n/magN, threshold);
506  }
507  else
508  {
510  << "Invalid normal for source face " << srcFacei
511  << " points " << UIndirectList<point>(srcPoints, src)
512  << " target face " << tgtFacei
513  << " points " << UIndirectList<point>(tgtPoints, tgt)
514  << endl;
515  }
516 
517  return false;
518 }
519 
520 
522 (
523  List<DynamicList<label>>& srcAddr,
524  List<DynamicList<scalar>>& srcWght,
525  List<DynamicList<point>>& srcCtr,
526  List<DynamicList<label>>& tgtAddr,
527  List<DynamicList<scalar>>& tgtWght
528 )
529 {
530  addProfiling(ami, "faceAreaWeightAMI::restartUncoveredSourceFace");
531 
532  // Note: exclude faces in srcNonOverlap_ for ACMI?
533 
534  label nBelowMinWeight = 0;
535  const scalar minWeight = 0.95;
536 
537  // List of tgt face neighbour faces
538  DynamicList<label> nbrFaces(10);
539 
540  // List of faces currently visited for srcFacei to avoid multiple hits
541  DynamicList<label> visitedFaces(10);
542 
543  const auto& srcPatch = this->srcPatch();
544 
545  forAll(srcWght, srcFacei)
546  {
547  const scalar s = sum(srcWght[srcFacei]);
548  const scalar t = s/srcMagSf_[srcFacei];
549 
550  if (t < minWeight)
551  {
552  ++nBelowMinWeight;
553 
554  const face& f = srcPatch[srcFacei];
555 
556  forAll(f, fpi)
557  {
558  const label tgtFacei =
559  findTargetFace(srcFacei, srcAddr[srcFacei], fpi);
560 
561  if (tgtFacei != -1)
562  {
563  nbrFaces.clear();
564  visitedFaces = srcAddr[srcFacei];
565 
566  (void)processSourceFace
567  (
568  srcFacei,
569  tgtFacei,
570 
571  nbrFaces,
572  visitedFaces,
573 
574  srcAddr,
575  srcWght,
576  srcCtr,
577  tgtAddr,
578  tgtWght
579  );
580  }
581  }
582  }
583  }
584 
585  if (debug && nBelowMinWeight)
586  {
588  << "Restarted search on " << nBelowMinWeight
589  << " faces since sum of weights < " << minWeight
590  << endl;
591  }
592 }
593 
594 
595 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
596 
598 (
599  const dictionary& dict,
600  const bool reverseTarget
601 )
602 :
603  advancingFrontAMI(dict, reverseTarget),
604  restartUncoveredSourceFace_
605  (
606  dict.getOrDefault("restartUncoveredSourceFace", true)
607  )
608 {}
609 
610 
612 (
613  const bool requireMatch,
614  const bool reverseTarget,
615  const scalar lowWeightCorrection,
617  const bool restartUncoveredSourceFace
618 )
619 :
621  (
622  requireMatch,
623  reverseTarget,
624  lowWeightCorrection,
625  triMode
626  ),
627  restartUncoveredSourceFace_(restartUncoveredSourceFace)
628 {}
629 
630 
632 :
633  advancingFrontAMI(ami),
634  restartUncoveredSourceFace_(ami.restartUncoveredSourceFace_)
635 {}
636 
637 
638 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
639 
641 (
642  const primitivePatch& srcPatch,
643  const primitivePatch& tgtPatch,
644  const autoPtr<searchableSurface>& surfPtr
645 )
646 {
647  if (upToDate_)
648  {
649  return false;
650  }
651 
652  addProfiling(ami, "faceAreaWeightAMI::calculate");
653 
654  advancingFrontAMI::calculate(srcPatch, tgtPatch, surfPtr);
655 
656  label srcFacei = 0;
657  label tgtFacei = 0;
658 
659  bool ok = initialiseWalk(srcFacei, tgtFacei);
660 
661  srcCentroids_.setSize(srcAddress_.size());
662 
663  const auto& src = this->srcPatch();
664  const auto& tgt = this->tgtPatch(); // might be the extended patch!
665 
666  // Temporary storage for addressing and weights
667  List<DynamicList<label>> srcAddr(src.size());
668  List<DynamicList<scalar>> srcWght(srcAddr.size());
669  List<DynamicList<point>> srcCtr(srcAddr.size());
670  List<DynamicList<label>> tgtAddr(tgt.size());
671  List<DynamicList<scalar>> tgtWght(tgtAddr.size());
672 
673  if (ok)
674  {
675  calcAddressing
676  (
677  srcAddr,
678  srcWght,
679  srcCtr,
680  tgtAddr,
681  tgtWght,
682  srcFacei,
683  tgtFacei
684  );
685 
686  if (debug && !srcNonOverlap_.empty())
687  {
688  Pout<< " AMI: " << srcNonOverlap_.size()
689  << " non-overlap faces identified"
690  << endl;
691  }
692 
693  // Check for badly covered faces
694  if (restartUncoveredSourceFace_) // && mustMatchFaces())
695  {
696  restartUncoveredSourceFace
697  (
698  srcAddr,
699  srcWght,
700  srcCtr,
701  tgtAddr,
702  tgtWght
703  );
704  }
705  }
706 
707  // Transfer data to persistent storage
708  forAll(srcAddr, i)
709  {
710  srcAddress_[i].transfer(srcAddr[i]);
711  srcWeights_[i].transfer(srcWght[i]);
712  srcCentroids_[i].transfer(srcCtr[i]);
713  }
714 
715  forAll(tgtAddr, i)
716  {
717  tgtAddress_[i].transfer(tgtAddr[i]);
718  tgtWeights_[i].transfer(tgtWght[i]);
719  }
720 
721  if (distributed())
722  {
723  const primitivePatch& srcPatch0 = this->srcPatch0();
724  const primitivePatch& tgtPatch0 = this->tgtPatch0();
725 
726  // Create global indexing for each original patch
727  globalIndex globalSrcFaces(srcPatch0.size());
728  globalIndex globalTgtFaces(tgtPatch0.size());
729 
730  for (labelList& addressing : srcAddress_)
731  {
732  for (label& addr : addressing)
733  {
734  addr = extendedTgtFaceIDs_[addr];
735  }
736  }
737 
738  for (labelList& addressing : tgtAddress_)
739  {
740  globalSrcFaces.inplaceToGlobal(addressing);
741  }
742 
743  // Send data back to originating procs. Note that contributions
744  // from different processors get added (ListOps::appendEqOp)
745 
747  (
749  List<labelPair>(),
750  tgtPatch0.size(),
751  extendedTgtMapPtr_->constructMap(),
752  false, // has flip
753  extendedTgtMapPtr_->subMap(),
754  false, // has flip
755  tgtAddress_,
756  labelList(),
758  flipOp() // flip operation
759  );
760 
762  (
764  List<labelPair>(),
765  tgtPatch0.size(),
766  extendedTgtMapPtr_->constructMap(),
767  false,
768  extendedTgtMapPtr_->subMap(),
769  false,
770  tgtWeights_,
771  scalarList(),
773  flipOp()
774  );
775 
776  // Note: using patch face areas calculated by the AMI method
777  extendedTgtMapPtr_->reverseDistribute(tgtPatch0.size(), tgtMagSf_);
778 
779  // Cache maps and reset addresses
780  List<Map<label>> cMapSrc;
781  srcMapPtr_.reset
782  (
783  new mapDistribute(globalSrcFaces, tgtAddress_, cMapSrc)
784  );
785 
786  List<Map<label>> cMapTgt;
787  tgtMapPtr_.reset
788  (
789  new mapDistribute(globalTgtFaces, srcAddress_, cMapTgt)
790  );
791  }
792 
793  // Convert the weights from areas to normalised values
794  normaliseWeights(requireMatch_, true);
795 
796  nonConformalCorrection();
797 
798  upToDate_ = true;
799 
800  return upToDate_;
801 }
802 
803 
805 {
807 
808  if (restartUncoveredSourceFace_)
809  {
810  os.writeEntry
811  (
812  "restartUncoveredSourceFace",
813  restartUncoveredSourceFace_
814  );
815  }
816 }
817 
818 
819 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
profiling.H
Foam::faceAreaWeightAMI::write
virtual void write(Ostream &os) const
Write.
Definition: faceAreaWeightAMI.C:804
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::faceAreaWeightAMI::restartUncoveredSourceFace
virtual void restartUncoveredSourceFace(List< DynamicList< label >> &srcAddr, List< DynamicList< scalar >> &srcWght, List< DynamicList< point >> &srcCtr, List< DynamicList< label >> &tgtAddr, List< DynamicList< scalar >> &tgtWght)
Attempt to re-evaluate source faces that have not been included.
Definition: faceAreaWeightAMI.C:522
Foam::faceAreaWeightAMI::processSourceFace
virtual bool processSourceFace(const label srcFacei, const label tgtStartFacei, DynamicList< label > &nbrFaces, DynamicList< label > &visitedFaces, List< DynamicList< label >> &srcAddr, List< DynamicList< scalar >> &srcWght, List< DynamicList< point >> &srcCtr, List< DynamicList< label >> &tgtAddr, List< DynamicList< scalar >> &tgtWght)
Determine overlap contributions for source face srcFacei.
Definition: faceAreaWeightAMI.C:183
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::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::OBJstream
OFstream that keeps track of vertices.
Definition: OBJstream.H:58
Foam::faceAreaIntersect::triangulationMode
triangulationMode
Definition: faceAreaIntersect.H:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< label >
Foam::bitSet::unset
bitSet & unset(const bitSet &other)
Definition: bitSetI.H:612
Foam::advancingFrontAMI::calculate
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
Definition: advancingFrontAMI.C:437
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::globalIndex::inplaceToGlobal
void inplaceToGlobal(labelList &labels) const
From local to global index (inplace)
Definition: globalIndexI.H:283
Foam::faceAreaIntersect::calc
void calc(const face &faceA, const face &faceB, const vector &n, scalar &area, vector &centroid) const
Definition: faceAreaIntersect.C:451
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::bitSet::test
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:520
Foam::OBJstream::write
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
Foam::DynamicList::clear
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
Foam::bitSet::count
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:499
Foam::faceAreaWeightAMI::setNextFaces
virtual bool setNextFaces(label &startSeedi, label &srcFacei, label &tgtFacei, const bitSet &mapFlag, labelList &seedFaces, const DynamicList< label > &visitedFaces, const bool errorOnNotFound=true) const
Set the source and target seed faces.
Definition: faceAreaWeightAMI.C:256
Foam::faceAreaWeightAMI::calcInterArea
virtual void calcInterArea(const label srcFacei, const label tgtFacei, scalar &area, vector &centroid) const
Area of intersection between source and target faces.
Definition: faceAreaWeightAMI.C:373
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::faceAreaWeightAMI
Face area weighted Arbitrary Mesh Interface (AMI) method.
Definition: faceAreaWeightAMI.H:53
Foam::faceAreaWeightAMI::faceAreaWeightAMI
faceAreaWeightAMI(const dictionary &dict, const bool reverseTarget=false)
Construct from dictionary.
Definition: faceAreaWeightAMI.C:598
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< vector >
Foam::mapDistributeBase::distribute
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &, const negateOp &negOp, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Distribute data. Note:schedule only used for.
Definition: mapDistributeBaseTemplates.C:122
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::faceAreaIntersect
Face intersection class.
Definition: faceAreaIntersect.H:59
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::faceAreaWeightAMI::calculate
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
Definition: faceAreaWeightAMI.C:641
Foam::AMIInterpolation::write
virtual void write(Ostream &os) const
Write.
Definition: AMIInterpolation.C:1243
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
dict
dictionary dict
Definition: searchingEngine.H:14
addProfiling
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
Definition: profilingTrigger.H:113
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::bitSet::find_next
label find_next(label pos) const
Locate the next bit set, starting one beyond the specified position.
Definition: bitSetI.H:400
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::faceAreaIntersect::overlaps
bool overlaps(const face &faceA, const face &faceB, const vector &n, const scalar threshold) const
Return area of intersection of faceA with faceB.
Definition: faceAreaIntersect.C:512
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::addAliasToRunTimeSelectionTable
addAliasToRunTimeSelectionTable(AMIInterpolation, faceAreaWeightAMI, dict, faceAreaWeightAMI, partialFaceAreaWeightAMI, 2012)
Foam::ListOps::appendEqOp
List helper to append y elements onto the end of x.
Definition: ListOps.H:570
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::UPstream::commsTypes::nonBlocking
f
labelList f(nPoints)
faceAreaWeightAMI.H
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::advancingFrontAMI
Base class for Arbitrary Mesh Interface (AMI) methods.
Definition: advancingFrontAMI.H:55
Foam::DynamicList::remove
T remove()
Remove and return the last element. Fatal on an empty list.
Definition: DynamicListI.H:704
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::PrimitivePatch::faceNormals
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
Definition: PrimitivePatch.C:445
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::PrimitivePatch::faceFaces
const labelListList & faceFaces() const
Return face-face addressing.
Definition: PrimitivePatch.C:249
Foam::faceAreaWeightAMI::calcAddressing
virtual void calcAddressing(List< DynamicList< label >> &srcAddress, List< DynamicList< scalar >> &srcWeights, List< DynamicList< point >> &srcCentroids, List< DynamicList< label >> &tgtAddress, List< DynamicList< scalar >> &tgtWeights, label srcFacei, label tgtFacei)
Definition: faceAreaWeightAMI.C:98
DebugVar
#define DebugVar(var)
Report a variable name and value.
Definition: messageStream.H:404
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::faceAreaIntersect::triangles
const DynamicList< triPoints > triangles() const
Const access to the triangulation.
Definition: faceAreaIntersectI.H:140
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::faceAreaWeightAMI::overlaps
virtual bool overlaps(const label srcFacei, const label tgtFacei, const scalar threshold) const
Return true if faces overlap.
Definition: faceAreaWeightAMI.C:456
OBJstream.H
Foam::triPointRef
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:71
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
Foam::flipOp
Functor to negate primitives. Dummy for most other types.
Definition: flipOp.H:53