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-2020 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 
43 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
44 
45 /*
46  if (debug)
47  {
48  static label nAMI = 0;
49 
50  // Write out triangulated surfaces as OBJ files
51  OBJstream srcTriObj("srcTris_" + Foam::name(nAMI) + ".obj");
52  const pointField& srcPts = src.points();
53  forAll(srcTris_, facei)
54  {
55  const DynamicList<face>& faces = srcTris_[facei];
56  for (const face& f : faces)
57  {
58  srcTriObj.write
59  (
60  triPointRef(srcPts[f[0]], srcPts[f[1]], srcPts[f[2]])
61  );
62  }
63  }
64 
65  OBJstream tgtTriObj("tgtTris_" + Foam::name(nAMI) + ".obj");
66  const pointField& tgtPts = tgt.points();
67  forAll(tgtTris_, facei)
68  {
69  const DynamicList<face>& faces = tgtTris_[facei];
70  for (const face& f : faces)
71  {
72  tgtTriObj.write
73  (
74  triPointRef(tgtPts[f[0]], tgtPts[f[1]], tgtPts[f[2]])
75  );
76  }
77  }
78 
79  ++nAMI;
80  }
81 */
82 
83 
85 (
86  List<DynamicList<label>>& srcAddr,
87  List<DynamicList<scalar>>& srcWght,
88  List<DynamicList<point>>& srcCtr,
89  List<DynamicList<label>>& tgtAddr,
90  List<DynamicList<scalar>>& tgtWght,
91  label srcFacei,
92  label tgtFacei
93 )
94 {
95  addProfiling(ami, "faceAreaWeightAMI::calcAddressing");
96 
97  // Construct weights and addressing
98  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
99 
100  label nFacesRemaining = srcAddr.size();
101 
102  // List of tgt face neighbour faces
103  DynamicList<label> nbrFaces(10);
104 
105  // List of faces currently visited for srcFacei to avoid multiple hits
106  DynamicList<label> visitedFaces(10);
107 
108  // List to keep track of tgt faces used to seed src faces
109  labelList seedFaces(nFacesRemaining, -1);
110  seedFaces[srcFacei] = tgtFacei;
111 
112  // List to keep track of whether src face can be mapped
113  bitSet mapFlag(nFacesRemaining, true);
114 
115  // Reset starting seed
116  label startSeedi = 0;
117 
118  bool continueWalk = true;
119  DynamicList<label> nonOverlapFaces;
120  do
121  {
122  nbrFaces.clear();
123  visitedFaces.clear();
124 
125  // Do advancing front starting from srcFacei,tgtFacei
126  bool faceProcessed = processSourceFace
127  (
128  srcFacei,
129  tgtFacei,
130 
131  nbrFaces,
132  visitedFaces,
133 
134  srcAddr,
135  srcWght,
136  srcCtr,
137  tgtAddr,
138  tgtWght
139  );
140 
141  mapFlag.unset(srcFacei);
142 
143  if (!faceProcessed)
144  {
145  nonOverlapFaces.append(srcFacei);
146  }
147 
148  // Choose new src face from current src face neighbour
149  continueWalk = setNextFaces
150  (
151  startSeedi,
152  srcFacei,
153  tgtFacei,
154  mapFlag,
155  seedFaces,
156  visitedFaces,
157  requireMatch_ && (lowWeightCorrection_ < 0)
158  // pass in nonOverlapFaces for failed tree search?
159  );
160  } while (continueWalk);
161 
162  srcNonOverlap_.transfer(nonOverlapFaces);
163 }
164 
165 
167 (
168  const label srcFacei,
169  const label tgtStartFacei,
170 
171  // list of tgt face neighbour faces
172  DynamicList<label>& nbrFaces,
173  // list of faces currently visited for srcFacei to avoid multiple hits
174  DynamicList<label>& visitedFaces,
175 
176  // temporary storage for addressing, weights and centroid
177  List<DynamicList<label>>& srcAddr,
178  List<DynamicList<scalar>>& srcWght,
179  List<DynamicList<point>>& srcCtr,
180  List<DynamicList<label>>& tgtAddr,
181  List<DynamicList<scalar>>& tgtWght
182 )
183 {
184  addProfiling(ami, "faceAreaWeightAMI::processSourceFace");
185 
186  if (tgtStartFacei == -1)
187  {
188  return false;
189  }
190 
191  const auto& tgtPatch = this->tgtPatch();
192 
193  // append initial target face and neighbours
194  nbrFaces.append(tgtStartFacei);
195  appendNbrFaces(tgtStartFacei, tgtPatch, visitedFaces, nbrFaces);
196 
197  bool faceProcessed = false;
198 
199  label maxNeighbourFaces = nbrFaces.size();
200 
201  do
202  {
203  // process new target face
204  label tgtFacei = nbrFaces.remove();
205  visitedFaces.append(tgtFacei);
206 
207  scalar interArea = 0;
208  vector interCentroid(Zero);
209  calcInterArea(srcFacei, tgtFacei, interArea, interCentroid);
210 
211  // store when intersection fractional area > tolerance
212  if (interArea/srcMagSf_[srcFacei] > faceAreaIntersect::tolerance())
213  {
214  srcAddr[srcFacei].append(tgtFacei);
215  srcWght[srcFacei].append(interArea);
216  srcCtr[srcFacei].append(interCentroid);
217 
218  tgtAddr[tgtFacei].append(srcFacei);
219  tgtWght[tgtFacei].append(interArea);
220 
221  appendNbrFaces(tgtFacei, tgtPatch, visitedFaces, nbrFaces);
222 
223  faceProcessed = true;
224 
225  maxNeighbourFaces = max(maxNeighbourFaces, nbrFaces.size());
226  }
227 
228  } while (nbrFaces.size() > 0);
229 
230  if (debug > 1)
231  {
232  DebugVar(maxNeighbourFaces);
233  }
234 
235  return faceProcessed;
236 }
237 
238 
240 (
241  label& startSeedi,
242  label& srcFacei,
243  label& tgtFacei,
244  const bitSet& mapFlag,
245  labelList& seedFaces,
246  const DynamicList<label>& visitedFaces,
247  const bool errorOnNotFound
248 ) const
249 {
250  addProfiling(ami, "faceAreaWeightAMI::setNextFaces");
251 
252  if (mapFlag.count() == 0)
253  {
254  // No more faces to map
255  return false;
256  }
257 
258  const labelList& srcNbrFaces = this->srcPatch().faceFaces()[srcFacei];
259 
260  // Initialise tgtFacei
261  tgtFacei = -1;
262 
263  // Set possible seeds for later use
264  bool valuesSet = false;
265  for (label faceS: srcNbrFaces)
266  {
267  if (mapFlag.test(faceS) && seedFaces[faceS] == -1)
268  {
269  for (label faceT : visitedFaces)
270  {
271  const scalar threshold =
272  srcMagSf_[faceS]*faceAreaIntersect::tolerance();
273 
274  // Store when intersection fractional area > threshold
275  if (overlaps(faceS, faceT, threshold))
276  {
277  seedFaces[faceS] = faceT;
278 
279  if (!valuesSet)
280  {
281  srcFacei = faceS;
282  tgtFacei = faceT;
283  valuesSet = true;
284  }
285  }
286  }
287  }
288  }
289 
290  if (valuesSet)
291  {
292  return true;
293  }
294 
295  // Set next src and tgt faces if not set above
296  // - try to use existing seed
297  label facei = startSeedi;
298  if (!mapFlag.test(startSeedi))
299  {
300  facei = mapFlag.find_next(facei);
301  }
302  const label startSeedi0 = facei;
303 
304  bool foundNextSeed = false;
305  while (facei != -1)
306  {
307  if (!foundNextSeed)
308  {
309  startSeedi = facei;
310  foundNextSeed = true;
311  }
312 
313  if (seedFaces[facei] != -1)
314  {
315  srcFacei = facei;
316  tgtFacei = seedFaces[facei];
317 
318  return true;
319  }
320 
321  facei = mapFlag.find_next(facei);
322  }
323 
324  // Perform new search to find match
325  if (debug)
326  {
327  Pout<< "Advancing front stalled: searching for new "
328  << "target face" << endl;
329  }
330 
331  facei = startSeedi0;
332  while (facei != -1)
333  {
334  srcFacei = facei;
335  tgtFacei = findTargetFace(srcFacei, visitedFaces);
336 
337  if (tgtFacei >= 0)
338  {
339  return true;
340  }
341 
342  facei = mapFlag.find_next(facei);
343  }
344 
345  if (errorOnNotFound)
346  {
348  << "Unable to set target face for source face " << srcFacei
349  << abort(FatalError);
350  }
351 
352  return false;
353 }
354 
355 
357 (
358  const label srcFacei,
359  const label tgtFacei,
360  scalar& area,
361  vector& centroid
362 ) const
363 {
364  addProfiling(ami, "faceAreaWeightAMI::interArea");
365 
366  // Quick reject if either face has zero area
367  if
368  (
369  (srcMagSf_[srcFacei] < ROOTVSMALL)
370  || (tgtMagSf_[tgtFacei] < ROOTVSMALL)
371  )
372  {
373  return;
374  }
375 
376  const auto& srcPatch = this->srcPatch();
377  const auto& tgtPatch = this->tgtPatch();
378 
379  const pointField& srcPoints = srcPatch.points();
380  const pointField& tgtPoints = tgtPatch.points();
381 
382  // Create intersection object
383  faceAreaIntersect inter
384  (
385  srcPoints,
386  tgtPoints,
387  srcTris_[srcFacei],
388  tgtTris_[tgtFacei],
389  reverseTarget_,
390  AMIInterpolation::cacheIntersections_
391  );
392 
393  // Crude resultant norm
394  vector n(-srcPatch.faceNormals()[srcFacei]);
395  if (reverseTarget_)
396  {
397  n -= tgtPatch.faceNormals()[tgtFacei];
398  }
399  else
400  {
401  n += tgtPatch.faceNormals()[tgtFacei];
402  }
403  scalar magN = mag(n);
404 
405  const face& src = srcPatch[srcFacei];
406  const face& tgt = tgtPatch[tgtFacei];
407 
408  if (magN > ROOTVSMALL)
409  {
410  inter.calc(src, tgt, n/magN, area, centroid);
411  }
412  else
413  {
415  << "Invalid normal for source face " << srcFacei
416  << " points " << UIndirectList<point>(srcPoints, src)
417  << " target face " << tgtFacei
418  << " points " << UIndirectList<point>(tgtPoints, tgt)
419  << endl;
420  }
421 
422  if (AMIInterpolation::cacheIntersections_ && debug)
423  {
424  static OBJstream tris("intersectionTris.obj");
425  const auto& triPts = inter.triangles();
426  for (const auto& tp : triPts)
427  {
428  tris.write(triPointRef(tp[0], tp[1], tp[2]), false);
429  }
430  }
431 
432  if ((debug > 1) && (area > 0))
433  {
434  writeIntersectionOBJ(area, src, tgt, srcPoints, tgtPoints);
435  }
436 }
437 
438 
440 (
441  const label srcFacei,
442  const label tgtFacei,
443  const scalar threshold
444 ) const
445 {
446  // Quick reject if either face has zero area
447  if
448  (
449  (srcMagSf_[srcFacei] < ROOTVSMALL)
450  || (tgtMagSf_[tgtFacei] < ROOTVSMALL)
451  )
452  {
453  return false;
454  }
455 
456  const auto& srcPatch = this->srcPatch();
457  const auto& tgtPatch = this->tgtPatch();
458 
459  const pointField& srcPoints = srcPatch.points();
460  const pointField& tgtPoints = tgtPatch.points();
461 
462  faceAreaIntersect inter
463  (
464  srcPoints,
465  tgtPoints,
466  srcTris_[srcFacei],
467  tgtTris_[tgtFacei],
468  reverseTarget_,
469  AMIInterpolation::cacheIntersections_
470  );
471 
472  // Crude resultant norm
473  vector n(-srcPatch.faceNormals()[srcFacei]);
474  if (reverseTarget_)
475  {
476  n -= tgtPatch.faceNormals()[tgtFacei];
477  }
478  else
479  {
480  n += tgtPatch.faceNormals()[tgtFacei];
481  }
482  scalar magN = mag(n);
483 
484  const face& src = srcPatch[srcFacei];
485  const face& tgt = tgtPatch[tgtFacei];
486 
487  if (magN > ROOTVSMALL)
488  {
489  return inter.overlaps(src, tgt, n/magN, threshold);
490  }
491  else
492  {
494  << "Invalid normal for source face " << srcFacei
495  << " points " << UIndirectList<point>(srcPoints, src)
496  << " target face " << tgtFacei
497  << " points " << UIndirectList<point>(tgtPoints, tgt)
498  << endl;
499  }
500 
501  return false;
502 }
503 
504 
506 (
507  List<DynamicList<label>>& srcAddr,
508  List<DynamicList<scalar>>& srcWght,
509  List<DynamicList<point>>& srcCtr,
510  List<DynamicList<label>>& tgtAddr,
511  List<DynamicList<scalar>>& tgtWght
512 )
513 {
514  addProfiling(ami, "faceAreaWeightAMI::restartUncoveredSourceFace");
515 
516  // Note: exclude faces in srcNonOverlap_ for ACMI?
517 
518  label nBelowMinWeight = 0;
519  const scalar minWeight = 0.95;
520 
521  // List of tgt face neighbour faces
522  DynamicList<label> nbrFaces(10);
523 
524  // List of faces currently visited for srcFacei to avoid multiple hits
525  DynamicList<label> visitedFaces(10);
526 
527  const auto& srcPatch = this->srcPatch();
528 
529  forAll(srcWght, srcFacei)
530  {
531  const scalar s = sum(srcWght[srcFacei]);
532  const scalar t = s/srcMagSf_[srcFacei];
533 
534  if (t < minWeight)
535  {
536  ++nBelowMinWeight;
537 
538  const face& f = srcPatch[srcFacei];
539 
540  forAll(f, fpi)
541  {
542  const label tgtFacei =
543  findTargetFace(srcFacei, srcAddr[srcFacei], fpi);
544 
545  if (tgtFacei != -1)
546  {
547  nbrFaces.clear();
548  visitedFaces = srcAddr[srcFacei];
549 
550  (void)processSourceFace
551  (
552  srcFacei,
553  tgtFacei,
554 
555  nbrFaces,
556  visitedFaces,
557 
558  srcAddr,
559  srcWght,
560  srcCtr,
561  tgtAddr,
562  tgtWght
563  );
564  }
565  }
566  }
567  }
568 
569  if (debug && nBelowMinWeight)
570  {
572  << "Restarted search on " << nBelowMinWeight
573  << " faces since sum of weights < " << minWeight
574  << endl;
575  }
576 }
577 
578 
579 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
580 
582 (
583  const dictionary& dict,
584  const bool reverseTarget
585 )
586 :
587  advancingFrontAMI(dict, reverseTarget),
588  restartUncoveredSourceFace_
589  (
590  dict.getOrDefault("restartUncoveredSourceFace", true)
591  )
592 {}
593 
594 
596 (
597  const bool requireMatch,
598  const bool reverseTarget,
599  const scalar lowWeightCorrection,
601  const bool restartUncoveredSourceFace
602 )
603 :
605  (
606  requireMatch,
607  reverseTarget,
608  lowWeightCorrection,
609  triMode
610  ),
611  restartUncoveredSourceFace_(restartUncoveredSourceFace)
612 {}
613 
614 
616 :
617  advancingFrontAMI(ami),
618  restartUncoveredSourceFace_(ami.restartUncoveredSourceFace_)
619 {}
620 
621 
622 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
623 
625 (
626  const primitivePatch& srcPatch,
627  const primitivePatch& tgtPatch,
628  const autoPtr<searchableSurface>& surfPtr
629 )
630 {
631  if (upToDate_)
632  {
633  return false;
634  }
635 
636  addProfiling(ami, "faceAreaWeightAMI::calculate");
637 
638  advancingFrontAMI::calculate(srcPatch, tgtPatch, surfPtr);
639 
640  label srcFacei = 0;
641  label tgtFacei = 0;
642 
643  bool ok = initialiseWalk(srcFacei, tgtFacei);
644 
645  srcCentroids_.setSize(srcAddress_.size());
646 
647  const auto& src = this->srcPatch();
648  const auto& tgt = this->tgtPatch(); // might be the extended patch!
649 
650  // Temporary storage for addressing and weights
651  List<DynamicList<label>> srcAddr(src.size());
652  List<DynamicList<scalar>> srcWght(srcAddr.size());
653  List<DynamicList<point>> srcCtr(srcAddr.size());
654  List<DynamicList<label>> tgtAddr(tgt.size());
655  List<DynamicList<scalar>> tgtWght(tgtAddr.size());
656 
657  if (ok)
658  {
659  calcAddressing
660  (
661  srcAddr,
662  srcWght,
663  srcCtr,
664  tgtAddr,
665  tgtWght,
666  srcFacei,
667  tgtFacei
668  );
669 
670  if (debug && !srcNonOverlap_.empty())
671  {
672  Pout<< " AMI: " << srcNonOverlap_.size()
673  << " non-overlap faces identified"
674  << endl;
675  }
676 
677  // Check for badly covered faces
678  if (restartUncoveredSourceFace_) // && requireMatch_???
679  {
680  restartUncoveredSourceFace
681  (
682  srcAddr,
683  srcWght,
684  srcCtr,
685  tgtAddr,
686  tgtWght
687  );
688  }
689  }
690 
691  // Transfer data to persistent storage
692  forAll(srcAddr, i)
693  {
694  srcAddress_[i].transfer(srcAddr[i]);
695  srcWeights_[i].transfer(srcWght[i]);
696  srcCentroids_[i].transfer(srcCtr[i]);
697  }
698 
699  forAll(tgtAddr, i)
700  {
701  tgtAddress_[i].transfer(tgtAddr[i]);
702  tgtWeights_[i].transfer(tgtWght[i]);
703  }
704 
705  if (distributed())
706  {
707  const primitivePatch& srcPatch0 = this->srcPatch0();
708  const primitivePatch& tgtPatch0 = this->tgtPatch0();
709 
710  // Create global indexing for each original patch
711  globalIndex globalSrcFaces(srcPatch0.size());
712  globalIndex globalTgtFaces(tgtPatch0.size());
713 
714  for (labelList& addressing : srcAddress_)
715  {
716  for (label& addr : addressing)
717  {
718  addr = extendedTgtFaceIDs_[addr];
719  }
720  }
721 
722  for (labelList& addressing : tgtAddress_)
723  {
724  globalSrcFaces.inplaceToGlobal(addressing);
725  }
726 
727  // Send data back to originating procs. Note that contributions
728  // from different processors get added (ListOps::appendEqOp)
729 
731  (
733  List<labelPair>(),
734  tgtPatch0.size(),
735  extendedTgtMapPtr_->constructMap(),
736  false, // has flip
737  extendedTgtMapPtr_->subMap(),
738  false, // has flip
739  tgtAddress_,
740  labelList(),
742  flipOp() // flip operation
743  );
744 
746  (
748  List<labelPair>(),
749  tgtPatch0.size(),
750  extendedTgtMapPtr_->constructMap(),
751  false,
752  extendedTgtMapPtr_->subMap(),
753  false,
754  tgtWeights_,
755  scalarList(),
757  flipOp()
758  );
759 
760  // Note: using patch face areas calculated by the AMI method
761  extendedTgtMapPtr_->reverseDistribute(tgtPatch0.size(), tgtMagSf_);
762 
763  // Cache maps and reset addresses
764  List<Map<label>> cMapSrc;
765  srcMapPtr_.reset
766  (
767  new mapDistribute(globalSrcFaces, tgtAddress_, cMapSrc)
768  );
769 
770  List<Map<label>> cMapTgt;
771  tgtMapPtr_.reset
772  (
773  new mapDistribute(globalTgtFaces, srcAddress_, cMapTgt)
774  );
775  }
776 
777  // Convert the weights from areas to normalised values
778  normaliseWeights(conformal(), true);
779 
780  upToDate_ = true;
781 
782  return upToDate_;
783 }
784 
785 
787 {
789 
790  if (restartUncoveredSourceFace_)
791  {
792  os.writeEntry
793  (
794  "restartUncoveredSourceFace",
795  restartUncoveredSourceFace_
796  );
797  }
798 }
799 
800 
801 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
profiling.H
Foam::faceAreaWeightAMI::write
virtual void write(Ostream &os) const
Write.
Definition: faceAreaWeightAMI.C:786
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::UPstream::commsTypes::nonBlocking
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:506
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:167
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:416
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::globalIndex::inplaceToGlobal
void inplaceToGlobal(labelList &labels) const
From local to global index (inplace)
Definition: globalIndexI.H:207
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::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:240
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:357
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:51
Foam::faceAreaWeightAMI::faceAreaWeightAMI
faceAreaWeightAMI(const dictionary &dict, const bool reverseTarget=false)
Construct from dictionary.
Definition: faceAreaWeightAMI.C:582
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:474
Foam::faceAreaIntersect
Face intersection class.
Definition: faceAreaIntersect.H:59
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:348
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:625
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:121
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
Foam::PrimitivePatch::points
const Field< point_type > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:305
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::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:381
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:653
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::PrimitivePatch::faceNormals
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
Definition: PrimitivePatch.C:425
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:229
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:85
DebugVar
#define DebugVar(var)
Report a variable name and value.
Definition: messageStream.H:381
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:303
Foam::faceAreaWeightAMI::overlaps
virtual bool overlaps(const label srcFacei, const label tgtFacei, const scalar threshold) const
Return true if faces overlap.
Definition: faceAreaWeightAMI.C:440
OBJstream.H
Foam::triPointRef
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:73
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:85
Foam::flipOp
Functor to negate primitives. Dummy for most other types.
Definition: flipOp.H:53