trackingInverseDistanceCellCellStencil.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) 2017-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify i
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "OBJstream.H"
31 #include "Time.H"
32 #include "fvMeshSubset.H"
33 #include "globalIndex.H"
34 #include "oversetFvPatch.H"
36 #include "syncTools.H"
37 #include "treeBoundBoxList.H"
38 #include "voxelMeshSearch.H"
39 #include "dynamicOversetFvMesh.H"
40 #include "waveMethod.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace cellCellStencils
47 {
48  defineTypeNameAndDebug(trackingInverseDistance, 0);
50 }
51 }
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
56 (
57  const fvMesh& mesh,
58  const vector& smallVec,
59 
60  const boundBox& bb,
61  labelVector& nDivs,
63 
64  const labelList& cellMap,
65  labelList& patchCellTypes
66 )
67 {
68  // Mark all voxels that overlap the bounding box of any patch
69 
70  static bool hasWarned = false;
71 
72  const fvBoundaryMesh& pbm = mesh.boundary();
73 
74  patchTypes = patchCellType::OTHER;
75 
76  // Mark wall boundaries
77  forAll(pbm, patchi)
78  {
79  const fvPatch& fvp = pbm[patchi];
80  const labelList& fc = fvp.faceCells();
81 
82  if (!fvPatch::constraintType(fvp.type()))
83  {
84  //Info<< "Marking cells on proper patch " << fvp.name()
85  // << " with type " << fvp.type() << endl;
86  const polyPatch& pp = fvp.patch();
87  forAll(pp, i)
88  {
89  // Mark in overall patch types
90  patchCellTypes[cellMap[fc[i]]] = patchCellType::PATCH;
91 
92  // Mark in voxel mesh
93  boundBox faceBb(pp.points(), pp[i], false);
94  faceBb.grow(smallVec);
95 
96  if (bb.overlaps(faceBb))
97  {
99  (
100  patchTypes,
101  bb,
102  nDivs,
103  faceBb,
105  );
106  }
107  }
108  }
109  }
110 
111  // Override with overset boundaries
112  forAll(pbm, patchi)
113  {
114  const fvPatch& fvp = pbm[patchi];
115  const labelList& fc = fvp.faceCells();
116 
117  if (isA<oversetFvPatch>(fvp))
118  {
119  //Info<< "Marking cells on overset patch " << fvp.name() << endl;
120  const polyPatch& pp = fvp.patch();
121  const vectorField::subField faceCentres(pp.faceCentres());
122  forAll(pp, i)
123  {
124  // Mark in overall patch types
125  patchCellTypes[cellMap[fc[i]]] = patchCellType::OVERSET;
126 
127  // Mark in voxel mesh
128  boundBox faceBb(pp.points(), pp[i], false);
129  faceBb.grow(smallVec);
130 
131  if (!bb.contains(faceCentres[i]))
132  {
133  if (!hasWarned)
134  {
135  hasWarned = true;
136  WarningInFunction << "Patch " << pp.name()
137  << " : face at " << faceCentres[i]
138  << " is not inside search bounding box of"
139  << " voxel mesh " << bb << endl
140  << " Is your 'searchBox' specification correct?"
141  << endl;
142  }
143  }
144  else
145  {
146  // Test for having voxels already marked as patch
147  // -> voxel mesh is too coarse
148  if
149  (
151  (
152  bb,
153  nDivs,
154  faceBb,
155  patchTypes,
156  static_cast<unsigned int>(patchCellType::PATCH)
157  )
158  )
159  {
160  // Determine number of voxels from number of cells
161  // in mesh
162  const labelVector& dim = mesh.geometricD();
163  forAll(dim, i)
164  {
165  if (dim[i] != -1)
166  {
167  nDivs[i] *= 2;
168  }
169  }
170 
171  Pout<< "Voxel mesh too coarse. Bounding box "
172  << faceBb
173  << " contains both non-overset and overset patches"
174  << ". Refining voxel mesh to " << nDivs << endl;
175 
176  return false;
177  }
178 
180  (
181  patchTypes,
182  bb,
183  nDivs,
184  faceBb,
185  patchCellType::OVERSET
186  );
187  }
188  }
189  }
190  }
191 
192  return true;
193 }
194 
195 
197 (
198  PstreamBuffers& pBufs,
199  const List<treeBoundBoxList>& patchBb,
200  const List<labelVector>& patchDivisions,
201  const PtrList<PackedList<2>>& patchParts,
202 
203  const label srcI,
204  const label tgtI,
205  labelList& allCellTypes
206 ) const
207 {
208  const treeBoundBoxList& srcPatchBbs = patchBb[srcI];
209  const treeBoundBoxList& tgtPatchBbs = patchBb[tgtI];
210  const labelList& tgtCellMap = meshParts_[tgtI].cellMap();
211 
212  // 1. do processor-local src-tgt patch overlap
213  {
214  const treeBoundBox& srcPatchBb = srcPatchBbs[Pstream::myProcNo()];
215  const treeBoundBox& tgtPatchBb = tgtPatchBbs[Pstream::myProcNo()];
216 
217  if (srcPatchBb.overlaps(tgtPatchBb))
218  {
219  const PackedList<2>& srcPatchTypes = patchParts[srcI];
220  const labelVector& srcDivs = patchDivisions[srcI];
221 
222  forAll(tgtCellMap, tgtCelli)
223  {
224  label celli = tgtCellMap[tgtCelli];
225  boundBox cBb(mesh_.cellBb(celli));
226  cBb.grow(smallVec_);
227 
228  if
229  (
231  (
232  srcPatchBb,
233  srcDivs,
234  cBb,
235  srcPatchTypes,
236  static_cast<unsigned int>(patchCellType::PATCH)
237  )
238  )
239  {
240  allCellTypes[celli] = HOLE;
241  }
242  }
243  }
244  }
245 
246 
247  // 2. Send over srcMesh bits that overlap tgt and do calculation
248  pBufs.clear();
249  for (const int procI : Pstream::allProcs())
250  {
251  if (procI != Pstream::myProcNo())
252  {
253  const treeBoundBox& srcPatchBb = srcPatchBbs[Pstream::myProcNo()];
254  const treeBoundBox& tgtPatchBb = tgtPatchBbs[procI];
255 
256  if (srcPatchBb.overlaps(tgtPatchBb))
257  {
258  // Send over complete patch voxel map. Tbd: could
259  // subset
260  UOPstream os(procI, pBufs);
261  os << srcPatchBb << patchDivisions[srcI] << patchParts[srcI];
262  }
263  }
264  }
265  pBufs.finishedSends();
266  for (const int procI : Pstream::allProcs())
267  {
268  if (procI != Pstream::myProcNo())
269  {
270  const treeBoundBox& srcPatchBb = srcPatchBbs[procI];
271  const treeBoundBox& tgtPatchBb = tgtPatchBbs[Pstream::myProcNo()];
272 
273  if (srcPatchBb.overlaps(tgtPatchBb))
274  {
275  UIPstream is(procI, pBufs);
276  const treeBoundBox receivedBb(is);
277  const labelVector srcDivs(is);
278  const PackedList<2> srcPatchTypes(is);
279 
280  // Verify validity
281  if (srcPatchBb != receivedBb)
282  {
284  << "proc:" << procI
285  << " srcPatchBb:" << srcPatchBb
286  << " receivedBb:" << receivedBb
287  << exit(FatalError);
288  }
289 
290  forAll(tgtCellMap, tgtCelli)
291  {
292  label celli = tgtCellMap[tgtCelli];
293  boundBox cBb(mesh_.cellBb(celli));
294  cBb.grow(smallVec_);
295 
296  if
297  (
299  (
300  srcPatchBb,
301  srcDivs,
302  cBb,
303  srcPatchTypes,
304  static_cast<unsigned int>(patchCellType::PATCH)
305  )
306  )
307  {
308  allCellTypes[celli] = HOLE;
309  }
310  }
311  }
312  }
313  }
314 }
315 
316 
318 (
319  PstreamBuffers& pBufs,
320  const List<treeBoundBoxList>& meshBb,
321  const PtrList<voxelMeshSearch>& meshSearches,
322  const labelList& allCellTypes,
323 
324  const label srcI,
325  const label tgtI,
326  labelListList& allStencil,
327  labelList& allDonor
328 ) const
329 {
330  const treeBoundBoxList& srcBbs = meshBb[srcI];
331  const treeBoundBoxList& tgtBbs = meshBb[tgtI];
332 
333  const fvMesh& srcMesh = meshParts_[srcI].subMesh();
334  const labelList& srcCellMap = meshParts_[srcI].cellMap();
335  const fvMesh& tgtMesh = meshParts_[tgtI].subMesh();
336  const pointField& tgtCc = tgtMesh.cellCentres();
337  const labelList& tgtCellMap = meshParts_[tgtI].cellMap();
338 
339  // 1. do processor-local src/tgt overlap
340  {
341  labelList tgtToSrcAddr;
342  waveMethod::calculate(tgtMesh, srcMesh, tgtToSrcAddr);
343 
344  forAll(tgtCellMap, tgtCelli)
345  {
346  const label srcCelli = tgtToSrcAddr[tgtCelli];
347 
348  if
349  (
350  srcCelli != -1
351  && allCellTypes[tgtCellMap[tgtCelli]] != HOLE
352  )
353  {
354  label celli = tgtCellMap[tgtCelli];
355 
356  // TBD: check for multiple donors. Maybe better one? For
357  // now check 'nearer' mesh
358  if (betterDonor(tgtI, allDonor[celli], srcI))
359  {
360  allStencil[celli].setSize(1);
361  allStencil[celli][0] =
362  globalCells_.toGlobal(srcCellMap[srcCelli]);
363  allDonor[celli] = srcI;
364  }
365  }
366  }
367  }
368 
369 
370  // 2. Send over tgtMesh bits that overlap src and do calculation on
371  // srcMesh.
372 
373 
374  // (remote) processors where the tgt overlaps my src
375  DynamicList<label> tgtOverlapProcs(Pstream::nProcs());
376  // (remote) processors where the src overlaps my tgt
377  DynamicList<label> srcOverlapProcs(Pstream::nProcs());
378  for (const int procI : Pstream::allProcs())
379  {
380  if (procI != Pstream::myProcNo())
381  {
382  if (tgtBbs[procI].overlaps(srcBbs[Pstream::myProcNo()]))
383  {
384  tgtOverlapProcs.append(procI);
385  }
386  if (srcBbs[procI].overlaps(tgtBbs[Pstream::myProcNo()]))
387  {
388  srcOverlapProcs.append(procI);
389  }
390  }
391  }
392 
393 
394 
395  // Indices of tgtcells to send over to each processor
396  List<DynamicList<label>> tgtSendCells(Pstream::nProcs());
397  forAll(srcOverlapProcs, i)
398  {
399  label procI = srcOverlapProcs[i];
400  tgtSendCells[procI].reserve(tgtMesh.nCells()/srcOverlapProcs.size());
401  }
402 
403 
404  forAll(tgtCellMap, tgtCelli)
405  {
406  label celli = tgtCellMap[tgtCelli];
407  if (srcOverlapProcs.size())
408  {
409  treeBoundBox subBb(mesh_.cellBb(celli));
410  subBb.grow(smallVec_);
411 
412  forAll(srcOverlapProcs, i)
413  {
414  label procI = srcOverlapProcs[i];
415  if (subBb.overlaps(srcBbs[procI]))
416  {
417  tgtSendCells[procI].append(tgtCelli);
418  }
419  }
420  }
421  }
422 
423  // Send target cell centres to overlapping processors
424  pBufs.clear();
425 
426  forAll(srcOverlapProcs, i)
427  {
428  label procI = srcOverlapProcs[i];
429  const labelList& cellIDs = tgtSendCells[procI];
430 
431  UOPstream os(procI, pBufs);
432  os << UIndirectList<point>(tgtCc, cellIDs);
433  }
434  pBufs.finishedSends();
435 
436  // Receive bits of target processors; find; send back
437  (void)srcMesh.tetBasePtIs();
438  forAll(tgtOverlapProcs, i)
439  {
440  label procI = tgtOverlapProcs[i];
441 
442  UIPstream is(procI, pBufs);
443  pointList samples(is);
444 
445  labelList donors(samples.size(), -1);
446  forAll(samples, sampleI)
447  {
448  const point& sample = samples[sampleI];
449  label srcCelli = srcMesh.findCell(sample, polyMesh::CELL_TETS);
450  if (srcCelli != -1)
451  {
452  donors[sampleI] = globalCells_.toGlobal(srcCellMap[srcCelli]);
453  }
454  }
455 
456  // Use same pStreamBuffers to send back.
457  UOPstream os(procI, pBufs);
458  os << donors;
459  }
460  pBufs.finishedSends();
461 
462  forAll(srcOverlapProcs, i)
463  {
464  label procI = srcOverlapProcs[i];
465  const labelList& cellIDs = tgtSendCells[procI];
466 
467  UIPstream is(procI, pBufs);
468  labelList donors(is);
469 
470  if (donors.size() != cellIDs.size())
471  {
472  FatalErrorInFunction<< "problem : cellIDs:" << cellIDs.size()
473  << " donors:" << donors.size() << abort(FatalError);
474  }
475 
476  forAll(donors, donorI)
477  {
478  label globalDonor = donors[donorI];
479 
480  if (globalDonor != -1)
481  {
482  label celli = tgtCellMap[cellIDs[donorI]];
483 
484  // TBD: check for multiple donors. Maybe better one?
485  if (betterDonor(tgtI, allDonor[celli], srcI))
486  {
487  allStencil[celli].setSize(1);
488  allStencil[celli][0] = globalDonor;
489  allDonor[celli] = srcI;
490  }
491  }
492  }
493  }
494 }
495 
496 
497 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
498 
499 Foam::cellCellStencils::trackingInverseDistance::trackingInverseDistance
500 (
501  const fvMesh& mesh,
502  const dictionary& dict,
503  const bool doUpdate
504 )
505 :
506  inverseDistance(mesh, dict, false),
507  globalCells_(mesh_.nCells())
508 {
509  // Initialise donor cell
511  globalDonor_ = -1;
512 
513  // Initialise mesh partitions
514  const labelIOList& zoneID = this->zoneID();
515  label nZones = gMax(zoneID)+1;
516 
517  labelList nCellsPerZone(nZones, Zero);
518  forAll(zoneID, celli)
519  {
520  nCellsPerZone[zoneID[celli]]++;
521  }
522  Pstream::listCombineReduce(nCellsPerZone, plusEqOp<label>());
523 
524  meshParts_.setSize(nZones);
525  forAll(meshParts_, zonei)
526  {
527  meshParts_.set
528  (
529  zonei,
530  new fvMeshSubset(mesh_, zonei, zoneID)
531  );
532 
533  // Trigger early evaluation of mesh dimension
534  // (in case there are locally zero cells in mesh)
535  (void)meshParts_[zonei].subMesh().nGeometricD();
536  }
537 
538 
539  // Print a bit
540  {
541  Info<< typeName << " : detected " << nZones
542  << " mesh regions" << endl;
543  Info<< incrIndent;
544  forAll(nCellsPerZone, zonei)
545  {
546  Info<< indent<< "zone:" << zonei
547  << " nCells:" << nCellsPerZone[zonei]
548  << endl;
549  }
550  Info<< decrIndent;
551  }
552 
553 
554  // Do geometry update
555  if (doUpdate)
556  {
558  }
559 }
560 
561 
562 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
565 {}
566 
567 
568 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
569 
571 {
572  DebugInfo<< FUNCTION_NAME << " : Start of analysis" << endl;
573 
574  scalar layerRelax(dict_.getOrDefault("layerRelax", 1.0));
575  const labelIOList& zoneID = this->zoneID();
576  label nZones = meshParts_.size();
577 
578  // Update stored mesh partitions for geometry changes
579  forAll(meshParts_, zonei)
580  {
581  pointField subPoints(mesh_.points(), meshParts_[zonei].pointMap());
582 
583  fvMesh& subMesh = meshParts_[zonei].subMesh();
584  subMesh.movePoints(subPoints);
585  }
586 
587  DebugInfo<< FUNCTION_NAME << " : Moved zone sub-meshes" << endl;
588 
589  // Calculate fast search structure for each zone
590  List<labelVector> searchBoxDivisions(dict_.lookup("searchBoxDivisions"));
591  if (searchBoxDivisions.size() != nZones)
592  {
593  FatalIOErrorInFunction(dict_) << "Number of searchBoxDivisions "
594  << searchBoxDivisions.size()
595  << " should equal the number of zones " << nZones
596  << exit(FatalIOError);
597  }
598 
599  const boundBox& allBb = mesh_.bounds();
600 
601  List<treeBoundBoxList> meshBb(nZones);
602 
603  // Determine zone meshes and bounding boxes
604  {
605  // Per processor, per zone the bounding box
606  List<treeBoundBoxList> procBb(Pstream::nProcs());
607  procBb[Pstream::myProcNo()].setSize(nZones);
608 
609  forAll(meshParts_, zonei)
610  {
611  const fvMesh& subMesh = meshParts_[zonei].subMesh();
612 
613  if (subMesh.nPoints())
614  {
615  procBb[Pstream::myProcNo()][zonei] =
616  treeBoundBox(subMesh.points());
617  procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
618  }
619  else
620  {
621  // No part of zone on this processor. Make up bb.
622  procBb[Pstream::myProcNo()][zonei] = treeBoundBox
623  (
624  allBb.min() - 2*allBb.span(),
625  allBb.min() - allBb.span()
626  );
627  procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
628  }
629  }
630 
631  Pstream::allGatherList(procBb);
632 
633  // Move local bounding boxes to per-mesh indexing
634  forAll(meshBb, zonei)
635  {
636  treeBoundBoxList& bbs = meshBb[zonei];
637  bbs.setSize(Pstream::nProcs());
638  forAll(procBb, proci)
639  {
640  bbs[proci] = procBb[proci][zonei];
641  }
642  }
643  }
644 
645  DebugInfo<< FUNCTION_NAME << " : Calculated bounding boxes" << endl;
646 
647 
648  // Determine patch bounding boxes. These are either global and provided
649  // by the user or processor-local as a copy of the mesh bounding box.
650 
651  List<treeBoundBoxList> patchBb(nZones);
652  List<labelVector> patchDivisions(searchBoxDivisions);
653  PtrList<PackedList<2>> patchParts(nZones);
654  labelList allPatchTypes(mesh_.nCells(), OTHER);
655 
656  {
657  treeBoundBox globalPatchBb;
658  if (dict_.readIfPresent("searchBox", globalPatchBb))
659  {
660  // All processors, all zones have the same bounding box
661  patchBb = treeBoundBoxList(Pstream::nProcs(), globalPatchBb);
662  }
663  else
664  {
665  // Use the meshBb (differing per zone, per processor)
666  patchBb = meshBb;
667  }
668  }
669 
670  forAll(patchParts, zonei)
671  {
672  while (true)
673  {
674  patchParts.set
675  (
676  zonei,
677  new PackedList<2>(cmptProduct(patchDivisions[zonei]))
678  );
679 
680  bool ok = markBoundaries
681  (
682  meshParts_[zonei].subMesh(),
683  smallVec_,
684 
685  patchBb[zonei][Pstream::myProcNo()],
686  patchDivisions[zonei],
687  patchParts[zonei],
688 
689  meshParts_[zonei].cellMap(),
690  allPatchTypes
691  );
692 
693  //if (returnReduceAnd(ok))
694  if (ok)
695  {
696  break;
697  }
698  }
699  }
700  DebugInfo<< FUNCTION_NAME << " : Calculated boundary voxel meshes" << endl;
701 
702 
703  PtrList<voxelMeshSearch> meshSearches(meshParts_.size());
704  forAll(meshParts_, zonei)
705  {
706  meshSearches.set
707  (
708  zonei,
709  new voxelMeshSearch
710  (
711  meshParts_[zonei].subMesh(),
712  patchBb[zonei][Pstream::myProcNo()],
713  patchDivisions[zonei],
714  true
715  )
716  );
717  }
718 
719  DebugInfo<< FUNCTION_NAME << " : Constructed cell voxel meshes" << endl;
720 
721  PtrList<fvMesh> voxelMeshes;
722  if (debug)
723  {
724  voxelMeshes.setSize(meshSearches.size());
725  forAll(meshSearches, zonei)
726  {
727  IOobject meshIO
728  (
729  word("voxelMesh" + Foam::name(zonei)),
730  mesh_.time().timeName(),
731  mesh_,
733  );
734 
735  Pout<< "Constructing voxel mesh " << meshIO.objectPath() << endl;
736  voxelMeshes.set(zonei, meshSearches[zonei].makeMesh(meshIO));
737  }
738  }
739 
740 
741  if (debug & 2)
742  {
743  forAll(patchParts, zonei)
744  {
745  const labelList voxelTypes(patchParts[zonei].values());
746  const fvMesh& vm = voxelMeshes[zonei];
747  tmp<volScalarField> tfld
748  (
749  createField<label>
750  (
751  vm,
752  "patchTypes",
753  voxelTypes
754  )
755  );
756  tfld().write();
757  }
758  }
759 
760  // Current best guess for cells
761  labelList allCellTypes(mesh_.nCells(), CALCULATED);
762  labelListList allStencil(mesh_.nCells());
763  // zoneID of donor
764  labelList allDonorID(mesh_.nCells(), -1);
765 
766  const globalIndex globalCells(mesh_.nCells());
767 
768  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
769 
770  DebugInfo<< FUNCTION_NAME << " : Allocated donor-cell structures" << endl;
771 
772  for (label srci = 0; srci < meshParts_.size()-1; srci++)
773  {
774  for (label tgti = srci+1; tgti < meshParts_.size(); tgti++)
775  {
776  markPatchesAsHoles
777  (
778  pBufs,
779 
780  patchBb,
781  patchDivisions,
782  patchParts,
783 
784  srci,
785  tgti,
786  allCellTypes
787  );
788  markPatchesAsHoles
789  (
790  pBufs,
791 
792  patchBb,
793  patchDivisions,
794  patchParts,
795 
796  tgti,
797  srci,
798  allCellTypes
799  );
800  }
801  }
802 
803  for (label srci = 0; srci < meshParts_.size()-1; srci++)
804  {
805  for (label tgti = srci+1; tgti < meshParts_.size(); tgti++)
806  {
807  markDonors
808  (
809  pBufs,
810  meshBb,
811  meshSearches,
812  allCellTypes, // to exclude hole donors
813 
814  tgti,
815  srci,
816  allStencil,
817  allDonorID
818  );
819  markDonors
820  (
821  pBufs,
822  meshBb,
823  meshSearches,
824  allCellTypes, // to exclude hole donors
825 
826  srci,
827  tgti,
828  allStencil,
829  allDonorID
830  );
831  }
832  }
833 
834  DebugInfo<< FUNCTION_NAME << " : Determined holes and donor-acceptors"
835  << endl;
836 
837  if ((debug & 2) && mesh_.time().writeTime())
838  {
839  tmp<volScalarField> tfld
840  (
841  createField(mesh_, "allCellTypes", allCellTypes)
842  );
843  tfld().write();
844 
845  tmp<volScalarField> tallDonorID
846  (
847  createField(mesh_, "allDonorID", allDonorID)
848  );
849  tallDonorID().write();
850  }
851 
852 
853  // Use the patch types and weights to decide what to do
854  forAll(allPatchTypes, celli)
855  {
856  if (allCellTypes[celli] != HOLE)
857  {
858  switch (allPatchTypes[celli])
859  {
860  case OVERSET:
861  {
862  // Require interpolation. See if possible.
863 
864  if (allStencil[celli].size())
865  {
866  allCellTypes[celli] = INTERPOLATED;
867  }
868  else
869  {
870  allCellTypes[celli] = HOLE;
871  }
872  }
873  }
874  }
875  }
876  DebugInfo<< FUNCTION_NAME << " : Removed bad donors" << endl;
877 
878  if ((debug & 2) && mesh_.time().writeTime())
879  {
880  tmp<volScalarField> tfld
881  (
882  createField(mesh_, "allCellTypes_patch", allCellTypes)
883  );
884  tfld().write();
885 
886  tmp<volScalarField> tfldOld
887  (
888  createField(mesh_, "allCellTypes_old", cellTypes_)
889  );
890  tfldOld().write();
891  }
892 
893  // Mark unreachable bits
894  findHoles(globalCells_, mesh_, zoneID, allStencil, allCellTypes);
895  DebugInfo<< FUNCTION_NAME << " : Flood-filled holes" << endl;
896 
897  if ((debug & 2) && mesh_.time().writeTime())
898  {
899  tmp<volScalarField> tfld
900  (
901  createField(mesh_, "allCellTypes_hole", allCellTypes)
902  );
903  tfld().write();
904 
905  labelList stencilSize(mesh_.nCells());
906  forAll(allStencil, celli)
907  {
908  stencilSize[celli] = allStencil[celli].size();
909  }
910  tmp<volScalarField> tfldsten
911  (
912  createField(mesh_, "allStencil_hole", stencilSize)
913  );
914  tfldsten().write();
915  }
916 
917  // Update allStencil with new fill HOLES
918  forAll(allCellTypes, celli)
919  {
920  if (allCellTypes[celli] == HOLE && allStencil[celli].size())
921  {
922  allStencil[celli].clear();
923  }
924  }
925 
926  // Add buffer interpolation layer(s) around holes
927  scalarField allWeight(mesh_.nCells(), Zero);
928 
929  labelListList compactStencil(allStencil);
930  List<Map<label>> compactStencilMap;
931  mapDistribute map(globalCells, compactStencil, compactStencilMap);
932 
933  scalarList compactCellVol(mesh_.V());
934  map.distribute(compactCellVol);
935 
936  walkFront
937  (
938  globalCells,
939  layerRelax,
940  allStencil,
941  allCellTypes,
942  allWeight,
943  compactCellVol,
944  compactStencil,
945  zoneID,
946  dict_.getOrDefault("holeLayers", 1),
947  dict_.getOrDefault("useLayer", -1)
948  );
949 
950  if ((debug & 2) && mesh_.time().writeTime())
951  {
952  tmp<volScalarField> tfld
953  (
954  createField(mesh_, "allCellTypes_front", allCellTypes)
955  );
956  tfld().write();
957  }
958 
959  // Check previous iteration cellTypes_ for any hole->calculated changes
960  // If so set the cell either to interpolated (if there are donors) or
961  // holes (if there are no donors). Note that any interpolated cell might
962  // still be overwritten by the flood filling
963  {
964  label nCalculated = 0;
965 
966  forAll(cellTypes_, celli)
967  {
968  if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
969  {
970  if (allStencil[celli].size() == 0)
971  {
972  // Reset to hole
973  allCellTypes[celli] = HOLE;
974  allStencil[celli].clear();
975  }
976  else
977  {
978  allCellTypes[celli] = INTERPOLATED;
979  nCalculated++;
980  }
981  }
982  }
983 
984  if (debug)
985  {
986  Pout<< "Detected " << nCalculated << " cells changing from hole"
987  << " to calculated. Changed to interpolated"
988  << endl;
989  }
990  }
991 
992 
993  // Convert cell-cell addressing to stencil in compact notation
994 
995  cellTypes_.transfer(allCellTypes);
996  cellStencil_.setSize(mesh_.nCells());
997  cellInterpolationWeights_.setSize(mesh_.nCells());
998  DynamicList<label> interpolationCells;
999 
1000  labelList compactCellTypes(cellTypes_);
1001  map.distribute(compactCellTypes);
1002 
1003  label nHoleDonors = 0;
1004  forAll(cellTypes_, celli)
1005  {
1006  if (cellTypes_[celli] == INTERPOLATED)
1007  {
1008  const labelList& slots = compactStencil[celli];
1009  if (slots.size())
1010  {
1011  if
1012  (
1013  compactCellTypes[slots[0]] == HOLE
1014  ||
1015  (
1016  !allowInterpolatedDonors_
1017  && compactCellTypes[slots[0]] == INTERPOLATED
1018  )
1019  )
1020  {
1021  cellTypes_[celli] = POROUS;
1022  cellStencil_[celli].clear();
1023  cellInterpolationWeights_[celli].clear();
1024  nHoleDonors++;
1025  }
1026  else
1027  {
1028  cellStencil_[celli].transfer(allStencil[celli]);
1029  cellInterpolationWeights_[celli].setSize(1);
1030  cellInterpolationWeights_[celli][0] = 1.0;
1031  interpolationCells.append(celli);
1032  }
1033  }
1034  }
1035  else
1036  {
1037  cellStencil_[celli].clear();
1038  cellInterpolationWeights_[celli].clear();
1039  }
1040  }
1041 
1042  reduce(nHoleDonors, sumOp<label>());
1043 
1044  DebugInfo<< FUNCTION_NAME << "nHole Donors : " << nHoleDonors << endl;
1045 
1046  interpolationCells_.transfer(interpolationCells);
1047 
1048  List<Map<label>> compactMap;
1049  cellInterpolationMap_.reset
1050  (
1051  new mapDistribute
1052  (
1053  globalCells,
1054  cellStencil_,
1055  compactMap
1056  )
1057  );
1058  cellInterpolationWeight_.transfer(allWeight);
1060  <
1062  oversetFvPatchField<scalar>
1063  >(cellInterpolationWeight_.boundaryFieldRef(), false);
1064 
1065 
1066  if ((debug & 2) && mesh_.time().writeTime())
1067  {
1068  // Dump mesh
1069  mesh_.time().write();
1070 
1071  // Dump stencil
1072  mkDir(mesh_.time().timePath());
1073  OBJstream str(mesh_.time().timePath()/"injectionStencil.obj");
1074  Pout<< typeName << " : dumping injectionStencil to "
1075  << str.name() << endl;
1076  pointField cc(mesh_.cellCentres());
1077  cellInterpolationMap().distribute(cc);
1078 
1079  forAll(cellStencil_, celli)
1080  {
1081  const labelList& slots = cellStencil_[celli];
1082  if (slots.size())
1083  {
1084  const point& accCc = mesh_.cellCentres()[celli];
1085  forAll(slots, i)
1086  {
1087  const point& donorCc = cc[slots[i]];
1088  str.writeLine(accCc, 0.1*accCc+0.9*donorCc);
1089  }
1090  }
1091  }
1092  }
1093 
1094  DebugInfo<< FUNCTION_NAME << " : Transferred donor to stencil" << endl;
1095 
1096 
1097  // Extend stencil to get inverse distance weighted neighbours
1098  createStencil(globalCells, allowHoleDonors_);
1099  DebugInfo<< FUNCTION_NAME << " : Extended stencil" << endl;
1100 
1101  // Optional: convert hole cells next to non-hole cells into
1102  // interpolate-from-neighbours (of cell type SPECIAL)
1103  if (allowHoleDonors_)
1104  {
1105  holeExtrapolationStencil(globalCells_);
1106  }
1107 
1108 
1109  if ((debug & 2) && mesh_.time().writeTime())
1110  {
1111  // Dump weight
1112  cellInterpolationWeight_.instance() = mesh_.time().timeName();
1113  cellInterpolationWeight_.write();
1114 
1115  // Dump stencil
1116  mkDir(mesh_.time().timePath());
1117  OBJstream str(mesh_.time().timePath()/"stencil.obj");
1118  Pout<< typeName << " : dumping to " << str.name() << endl;
1119  pointField cc(mesh_.cellCentres());
1120  cellInterpolationMap().distribute(cc);
1121 
1122  forAll(cellStencil_, celli)
1123  {
1124  const labelList& slots = cellStencil_[celli];
1125  if (slots.size())
1126  {
1127  const point& accCc = mesh_.cellCentres()[celli];
1128  forAll(slots, i)
1129  {
1130  const point& donorCc = cc[slots[i]];
1131  str.writeLine(accCc, 0.1*accCc+0.9*donorCc);
1132  }
1133  }
1134  }
1135  }
1136 
1137 
1138  // Print some stats
1139  Info<< this->info();
1140 
1141  DebugInfo<< FUNCTION_NAME << " : Finished analysis" << endl;
1142 
1143  // Tbd: detect if anything changed. Most likely it did!
1144  return true;
1145 }
1146 
1147 
1148 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
const polyBoundaryMesh & pbm
Inverse-distance-weighted interpolation stencil.
dictionary dict
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:597
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool contains(const point &pt) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:455
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
bool overlaps(const boundBox &bb) const
Overlaps with other bounding box, sphere etc?
Definition: boundBoxI.H:439
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:1176
wordList patchTypes(nPatches)
static void calculate(const polyMesh &src, const polyMesh &tgt, labelList &srcToTgtAddr)
Calculate addressing.
Definition: waveMethod.C:38
scalarField samples(nIntervals, Zero)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:918
defineTypeNameAndDebug(cellVolumeWeight, 0)
PtrList< fvMeshSubset > meshParts_
Subset according to zone.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
Definition: treeBoundBox.H:83
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
PDRpatchDef PATCH
Definition: PDRpatchDef.H:120
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Calculation of interpolation stencils.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:107
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static bool overlaps(const boundBox &bb, const labelVector &nDivs, const boundBox &subBb, const Container &elems, const Type val, const bool isNot=false)
Check if any voxel inside bounding box is set to val or.
void markPatchesAsHoles(PstreamBuffers &pBufs, const List< treeBoundBoxList > &patchBb, const List< labelVector > &patchDivisions, const PtrList< PackedList< 2 >> &patchParts, const label srcI, const label tgtI, labelList &allCellTypes) const
Mark all cells overlapping (a voxel covered by) a src patch.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const fvMesh & mesh_
Reference to the mesh.
#define DebugInfo
Report an information message using Foam::Info.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
static void correctBoundaryConditions(typename GeoField::Boundary &bfld, const bool typeOnly)
Correct boundary conditions of certain type (typeOnly = true)
virtual bool update()
Update stencils. Return false if nothing changed.
addToRunTimeSelectionTable(cellCellStencil, cellVolumeWeight, mesh)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
Buffers for inter-processor communications streams (UOPstream, UIPstream).
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: fvPatch.C:74
void clear()
Clear all send/recv buffers and reset states.
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
void markDonors(PstreamBuffers &pBufs, const List< treeBoundBoxList > &meshBb, const PtrList< voxelMeshSearch > &meshSearches, const labelList &allCellTypes, const label srcI, const label tgtI, labelListList &allStencil, labelList &allDonor) const
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
label nCells() const noexcept
Number of mesh cells.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
static void fill(Container &elems, const boundBox &bb, const labelVector &nDivs, const boundBox &subBb, const Type val)
Fill voxels indicated by bounding box.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
"nonBlocking" : (MPI_Isend, MPI_Irecv)
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
messageStream Info
Information stream (stdout output on master, null elsewhere)
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
Vector< label > labelVector
Vector of labels.
Definition: labelVector.H:47
List< point > pointList
List of point.
Definition: pointList.H:32
List< label > labelList
A List of labels.
Definition: List.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition: fvPatch.H:202
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:439
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
labelList cellIDs
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
Definition: boundBoxI.H:367
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
static bool markBoundaries(const fvMesh &mesh, const vector &smallVec, const boundBox &bb, labelVector &nDivs, PackedList< 2 > &patchTypes, const labelList &cellMap, labelList &patchCellTypes)
Mark voxels of patchTypes with type of patch face.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static tmp< volScalarField > createField(const fvMesh &mesh, const scalar val)
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127