cellCellStencil.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 
28 #include "cellCellStencil.H"
30 #include "volFields.H"
31 #include "syncTools.H"
32 #include "globalIndex.H"
33 #include "oversetFvPatchFields.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(cellCellStencil, 0);
40  defineRunTimeSelectionTable(cellCellStencil, mesh);
41 }
42 
43 const Foam::Enum
44 <
46 >
48 ({
49  { cellType::CALCULATED, "calculated" },
50  { cellType::INTERPOLATED, "interpolated" },
51  { cellType::HOLE, "hole" },
52  { cellType::SPECIAL, "special" },
53  { cellType::POROUS, "porous" },
54 });
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
59 Foam::cellCellStencil::cellCellStencil(const fvMesh& mesh)
60 :
61  mesh_(mesh),
62  nonInterpolatedFields_({"zoneID"})
63 {}
64 
65 
67 (
68  const fvMesh& mesh,
69  const dictionary& dict,
70  const bool update
71 )
72 {
73  DebugInFunction << "Constructing cellCellStencil" << endl;
74 
75  const word stencilType(dict.get<word>("method"));
76 
77  auto* ctorPtr = meshConstructorTable(stencilType);
78 
79  if (!ctorPtr)
80  {
82  (
83  dict,
84  "cellCellStencil",
85  stencilType,
86  *meshConstructorTablePtr_
87  ) << exit(FatalIOError);
88  }
89 
90  return autoPtr<cellCellStencil>(ctorPtr(mesh, dict, update));
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
95 
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 {
104  if (name.ends_with("_0"))
105  {
106  return baseName(name.substr(0, name.size() - 2));
107  }
108 
109  return name;
110 }
111 
112 
114 {
115  // Protect local fields from interpolation
116  nonInterpolatedFields_.insert("cellInterpolationWeight");
117  nonInterpolatedFields_.insert("cellTypes");
118  nonInterpolatedFields_.insert("maxMagWeight");
119 
120  // For convenience also suppress frequently used displacement field
121  {
122  nonInterpolatedFields_.insert("cellDisplacement");
123  nonInterpolatedFields_.insert("grad(cellDisplacement)");
124  const word w("snGradCorr(cellDisplacement)");
125  const word d("((viscosity*faceDiffusivity)*magSf)");
126  nonInterpolatedFields_.insert("surfaceIntegrate(("+d+"*"+w+"))");
127  }
128 
129  // For convenience also suppress frequently used displacement field
130  {
131  nonInterpolatedFields_.insert("cellMotionU");
132  nonInterpolatedFields_.insert("grad(cellMotionU)");
133  const word w("snGradCorr(cellMotionU)");
134  const word d("((viscosity*faceDiffusivity)*magSf)");
135  nonInterpolatedFields_.insert("surfaceIntegrate(("+d+"*"+w+"))");
136  }
137 }
138 
139 
141 {
142  labelIOList* zoneIDPtr = mesh.getObjectPtr<labelIOList>("zoneID");
143 
144  if (!zoneIDPtr)
145  {
146  zoneIDPtr = new labelIOList
147  (
148  IOobject
149  (
150  "zoneID",
153  mesh,
157  ),
158  mesh.nCells()
159  );
160 
161  zoneIDPtr->store();
162 
163  auto& zoneID = *zoneIDPtr;
164 
165  volScalarField volZoneID
166  (
167  IOobject
168  (
169  "zoneID",
170  mesh.time().findInstance(mesh.dbDir(), "zoneID"),
171  mesh,
175  ),
176  mesh
177  );
178  forAll(volZoneID, celli)
179  {
180  zoneID[celli] = label(volZoneID[celli]);
181  }
182  }
183 
184  return *zoneIDPtr;
185 }
186 
187 
189 (
190  const label size,
191  const labelUList& lst
192 )
193 {
194  labelList count(size, Zero);
195  forAll(lst, i)
196  {
197  count[lst[i]]++;
198  }
200  return count;
201 }
202 
205 {
206  return nonInterpolatedFields_;
207 }
208 
211 {
212  return nonInterpolatedFields_;
213 }
214 
215 
216 bool Foam::cellCellStencil::localStencil(const labelUList& slots) const
217 {
218  forAll(slots, i)
219  {
220  if (slots[i] >= mesh_.nCells())
221  {
222  return false;
223  }
224  }
225  return true;
226 }
227 
228 
230 (
231  const globalIndex& gi,
232  const polyMesh& mesh,
233  const boolList& isValidCell,
234  const labelList& selectedCells,
235  labelListList& cellCells,
236  pointListList& cellCellCentres
237 )
238 {
239  // For selected cells determine the face neighbours (in global numbering)
240 
241  const pointField& cellCentres = mesh.cellCentres();
242  const labelList& faceOwner = mesh.faceOwner();
243  const labelList& faceNeighbour = mesh.faceNeighbour();
244  const cellList& cells = mesh.cells();
245 
246 
247  // 1. Determine global cell number on other side of coupled patches
248 
249  labelList globalCellIDs(identity(gi.localSize(), gi.localStart()));
250 
251  labelList nbrGlobalCellIDs;
253  (
254  mesh,
255  globalCellIDs,
256  nbrGlobalCellIDs
257  );
258  pointField nbrCellCentres;
260  (
261  mesh,
262  cellCentres,
263  nbrCellCentres
264  );
265 
266  boolList nbrIsValidCell;
268  (
269  mesh,
270  isValidCell,
271  nbrIsValidCell
272  );
273 
274 
275  // 2. Collect cell and all its neighbours
276 
277  cellCells.setSize(mesh.nCells());
278  cellCellCentres.setSize(cellCells.size());
279 
280  forAll(selectedCells, i)
281  {
282  label celli = selectedCells[i];
283 
284  const cell& cFaces = cells[celli];
285  labelList& stencil = cellCells[celli];
286  pointList& stencilPoints = cellCellCentres[celli];
287  stencil.setSize(cFaces.size()+1);
288  stencilPoints.setSize(stencil.size());
289  label compacti = 0;
290 
291  // First entry is cell itself
292  if (isValidCell[celli])
293  {
294  stencil[compacti] = globalCellIDs[celli];
295  stencilPoints[compacti++] = cellCentres[celli];
296  }
297 
298  // Other entries are cell neighbours
299  forAll(cFaces, i)
300  {
301  label facei = cFaces[i];
302  label bFacei = facei-mesh.nInternalFaces();
303  label own = faceOwner[facei];
304  label nbrCelli;
305  point nbrCc;
306  bool isValid = false;
307  if (bFacei >= 0)
308  {
309  nbrCelli = nbrGlobalCellIDs[bFacei];
310  nbrCc = nbrCellCentres[bFacei];
311  isValid = nbrIsValidCell[bFacei];
312  }
313  else
314  {
315  if (own != celli)
316  {
317  nbrCelli = gi.toGlobal(own);
318  nbrCc = cellCentres[own];
319  isValid = isValidCell[own];
320  }
321  else
322  {
323  label nei = faceNeighbour[facei];
324  nbrCelli = gi.toGlobal(nei);
325  nbrCc = cellCentres[nei];
326  isValid = isValidCell[nei];
327  }
328  }
329 
330  if (isValid)
331  {
332  SubList<label> current(stencil, compacti);
333  if (!current.found(nbrCelli))
334  {
335  stencil[compacti] = nbrCelli;
336  stencilPoints[compacti++] = nbrCc;
337  }
338  }
339  }
340  stencil.setSize(compacti);
341  stencilPoints.setSize(compacti);
342  }
343 }
344 
346 (
347  const label celli,
348  const scalar wantedFraction,
349  bitSet& isFront,
350  scalarField& fraction
351 ) const
352 {
353  const cell& cFaces = mesh_.cells()[celli];
354  forAll(cFaces, i)
355  {
356  const label nbrFacei = cFaces[i];
357  if (fraction[nbrFacei] < wantedFraction)
358  {
359  fraction[nbrFacei] = wantedFraction;
360  isFront.set(nbrFacei);
361  }
362  }
363 }
364 
365 
367 (
368  const labelList& allCellTypes,
369  bitSet& isFront
370 ) const
371 {
372  const labelList& own = mesh_.faceOwner();
373  const labelList& nei = mesh_.faceNeighbour();
374 
375  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
376  {
377  const label ownType = allCellTypes[own[facei]];
378  const label neiType = allCellTypes[nei[facei]];
379  if
380  (
381  (ownType == HOLE && neiType != HOLE)
382  || (ownType != HOLE && neiType == HOLE)
383  )
384  {
385  isFront.set(facei);
386  }
387  }
388 
389  labelList nbrCellTypes;
390  syncTools::swapBoundaryCellList(mesh_, allCellTypes, nbrCellTypes);
391 
392  for
393  (
394  label facei = mesh_.nInternalFaces();
395  facei < mesh_.nFaces();
396  facei++
397  )
398  {
399  const label ownType = allCellTypes[own[facei]];
400  const label neiType = nbrCellTypes[facei-mesh_.nInternalFaces()];
401 
402  if
403  (
404  (ownType == HOLE && neiType != HOLE)
405  || (ownType != HOLE && neiType == HOLE)
406  )
407  {
408  isFront.set(facei);
409  }
410  }
411 }
412 
413 
415 (
416  const labelList& allCellTypes,
417  bitSet& isFront
418 ) const
419 {
420  const fvBoundaryMesh& fvm = mesh_.boundary();
421  // 'overset' patches
422  forAll(fvm, patchi)
423  {
424  if (isA<oversetFvPatch>(fvm[patchi]))
425  {
426  const labelList& fc = fvm[patchi].faceCells();
427  forAll(fc, i)
428  {
429  const label celli = fc[i];
430  if (allCellTypes[celli] == INTERPOLATED)
431  {
432  // Note that acceptors might have been marked hole if
433  // there are no donors in which case we do not want to
434  // walk this out. This is an extreme situation.
435  isFront.set(fvm[patchi].start()+i);
436  }
437  }
438  }
439  }
440 }
441 
442 
444 (
445  const globalIndex& globalCells,
446  const scalar layerRelax,
447  const labelListList& allStencil,
448  labelList& allCellTypes,
449  scalarField& allWeight,
450  const scalarList& compactCellVol,
451  const labelListList& compactStencil,
452  const labelList& zoneID,
453  const label holeLayers,
454  const label useLayer
455 ) const
456 {
457  if (useLayer > holeLayers)
458  {
459  FatalErrorInFunction<< "useLayer: " << useLayer
460  << " is larger than : " << holeLayers
461  << abort(FatalError);
462  }
463 
464  if (useLayer == 0)
465  {
466  FatalErrorInFunction<< "useLayer: " << useLayer
467  << " can not be zero."
468  << abort(FatalError);
469  }
470 
471  // Current front
472  bitSet isFront(mesh_.nFaces());
473 
474  // List of cellTYpes
475  DynamicList<labelList> dataSet;
476  // List of cellTypes for increasing layers
477  DynamicList<scalarField> dAlllWeight;
478  // List of average volumen ration interpolated/donor
479  DynamicList<scalar> daverageVolRatio;
480 
481  // Counting holes of different layers
482  DynamicList<label> dHoles;
483 
484  const scalarField& V = mesh_.V();
485 
486  // Current layer
487  labelField nLayer(allCellTypes.size(), 0);
488 
489  for (label currLayer = 1; currLayer < holeLayers+1; ++currLayer)
490  {
491  // Set up original front
492  isFront = false;
493 
494  setUpFront(allCellTypes, isFront);
495 
496  labelList allCellTypesWork(allCellTypes);
497 
498  bitSet isFrontWork(isFront);
499  label nCurrLayer = currLayer;
500 
501  while (nCurrLayer > 1 && returnReduce(isFrontWork.any(), orOp<bool>()))
502  {
503  bitSet newIsFront(mesh_.nFaces());
504  forAll(isFrontWork, facei)
505  {
506  if (isFrontWork.test(facei))
507  {
508  const label own = mesh_.faceOwner()[facei];
509 
510  if (allCellTypesWork[own] != HOLE)
511  {
512  allCellTypesWork[own] = HOLE;
513  newIsFront.set(mesh_.cells()[own]);
514  }
515 
516  if (mesh_.isInternalFace(facei))
517  {
518  const label nei = mesh_.faceNeighbour()[facei];
519  if (allCellTypesWork[nei] != HOLE)
520  {
521  allCellTypesWork[nei] = HOLE;
522  newIsFront.set(mesh_.cells()[nei]);
523  }
524  }
525  }
526  }
527  syncTools::syncFaceList(mesh_, newIsFront, orEqOp<unsigned int>());
528 
529  isFrontWork.transfer(newIsFront);
530 
531  nCurrLayer--;
532  }
533 
534 
535  if ((debug & 2) && mesh_.time().writeTime())
536  {
537  tmp<volScalarField> tfld
538  (
540  (
541  mesh_,
542  "allCellTypesWork_Holes" + name(currLayer),
543  allCellTypesWork
544  )
545  );
546  tfld().write();
547  }
548 
549  if (currLayer == 1)
550  {
551  setUpFrontOnOversetPatch(allCellTypes, isFront);
552  }
553 
554  if (currLayer > 1)
555  {
556  isFront = false;
557  setUpFrontOnOversetPatch(allCellTypesWork, isFront);
558  setUpFront(allCellTypesWork, isFront);
559  }
560 
561  // Current interpolation fraction
562  scalarField fraction(mesh_.nFaces(), Zero);
563 
564  // Ratio between Inter/donor
565  scalarField volRatio(allCellTypes.size(), Zero);
566 
567  forAll(isFront, facei)
568  {
569  if (isFront.test(facei))
570  {
571  fraction[facei] = 1.0;
572  }
573  }
574 
575  scalarField allWeightWork(allCellTypes.size(), Zero);
576  bitSet nHoles(allCellTypes.size());
577 
578  while (returnReduce(isFront.any(), orOp<bool>()))
579  {
580  // Interpolate cells on front
581  bitSet newIsFront(mesh_.nFaces());
582  scalarField newFraction(fraction);
583 
584  forAll(isFront, facei)
585  {
586  if (isFront.test(facei))
587  {
588  const label own = mesh_.faceOwner()[facei];
589 
590  if (allCellTypesWork[own] != HOLE)
591  {
592  if (allWeightWork[own] < fraction[facei])
593  {
594  // Cell wants to become interpolated (if sufficient
595  // stencil, otherwise becomes hole)
596  if (allStencil[own].size())
597  {
598  nLayer[own] = currLayer;
599 
600  allWeightWork[own] = fraction[facei];
601  allCellTypesWork[own] = INTERPOLATED;
602 
603  const label donorId = compactStencil[own][0];
604 
605  volRatio[own] = V[own]/compactCellVol[donorId];
606 
607  seedCell
608  (
609  own,
610  fraction[facei]-layerRelax,
611  newIsFront,
612  newFraction
613  );
614 
615  nHoles[own] = true;
616  }
617  else
618  {
619  allWeightWork[own] = 0.0;
620  allCellTypesWork[own] = HOLE;
621  // Add faces of cell as new front
622  seedCell
623  (
624  own,
625  1.0,
626  newIsFront,
627  newFraction
628  );
629  }
630  }
631  }
632  if (mesh_.isInternalFace(facei))
633  {
634  label nei = mesh_.faceNeighbour()[facei];
635  if (allCellTypesWork[nei] != HOLE)
636  {
637  if (allWeightWork[nei] < fraction[facei])
638  {
639  if (allStencil[nei].size())
640  {
641  nLayer[nei] = currLayer;
642 
643  allWeightWork[nei] = fraction[facei];
644  allCellTypesWork[nei] = INTERPOLATED;
645 
646  const label donorId =
647  compactStencil[nei][0];
648 
649  volRatio[nei] =
650  V[nei]/compactCellVol[donorId];
651 
652  seedCell
653  (
654  nei,
655  fraction[facei]-layerRelax,
656  newIsFront,
657  newFraction
658  );
659  }
660  else
661  {
662  allWeightWork[nei] = 0.0;
663  allCellTypesWork[nei] = HOLE;
664  nHoles[nei] = true;
665  seedCell
666  (
667  nei,
668  1.0,
669  newIsFront,
670  newFraction
671  );
672  }
673  }
674  }
675  }
676  }
677  }
678 
679  syncTools::syncFaceList(mesh_, newIsFront, orEqOp<unsigned int>());
680  syncTools::syncFaceList(mesh_, newFraction, maxEqOp<scalar>());
681 
682  isFront.transfer(newIsFront);
683  fraction.transfer(newFraction);
684  }
685 
686  if ((debug & 2) && mesh_.time().writeTime())
687  {
688  tmp<volScalarField> tfld
689  (
691  (
692  mesh_,
693  "allCellTypesWork_Layers" + name(currLayer),
694  allCellTypesWork
695  )
696  );
697  tfld().write();
698  }
699 
700  dAlllWeight.append(allWeightWork);
701  dataSet.append(allCellTypesWork);
702 
703  // Counting interpolated cells
704  label count = 1;
705  forAll (volRatio, i)
706  {
707  if (volRatio[i] > 0)
708  {
709  count++;
710  }
711  }
712  label nCount(count);
713  reduce(nCount, sumOp<label>());
714 
715  const scalar aveVol = mag(gSum(volRatio));
716 
717  daverageVolRatio.append(aveVol/nCount);
718 
719  // Check holes number. A sudden increase occurs when the walk leaks
720  // out from the obstacle
721  label nTotalHoles(nHoles.count());
722  reduce(nTotalHoles, sumOp<label>());
723  dHoles.append(nTotalHoles);
724 
725  // Check the increase between this layer and the last one
726  // if over 50% breaks the layer loop.
727  if
728  (
729  (currLayer > 1)
730  & (nTotalHoles > 2.0*dHoles[currLayer - 1])
731  )
732  {
733  break;
734  }
735  }
736 
737 
738  if ((debug & 2) && mesh_.time().writeTime())
739  {
740  tmp<volScalarField> tfld
741  (
742  createField(mesh_, "walkFront_layers", nLayer)
743  );
744  tfld().write();
745  }
746 
747  // Try to find the best averageVolRatio the further from the initial set
748  // As this one is next to HOLES
749  scalarList averageVolRatio;
750  averageVolRatio.transfer(daverageVolRatio);
751  label minVolId = findMin(averageVolRatio);
752 
753  if (useLayer > -1)
754  {
755  Info<< nl << " Number of layers : " << averageVolRatio.size() << nl
756  << " Average volumetric ratio : " << averageVolRatio << nl
757  << " Number of holes cells : " << dHoles << nl
758  << " Using layer : " << useLayer << nl << endl;
759  }
760 
761  if (useLayer != -1)
762  {
763  minVolId = useLayer - 1;
764  }
765 
766  allCellTypes.transfer(dataSet[minVolId]);
767  allWeight.transfer(dAlllWeight[minVolId]);
768 }
769 
770 
771 // * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * * * //
772 
773 template<>
774 Foam::Ostream& Foam::operator<<
775 (
776  Ostream& os,
777  const InfoProxy<cellCellStencil>& iproxy
778 )
779 {
780  const auto& e = *iproxy;
781 
782  const labelUList& cellTypes = e.cellTypes();
783  const labelUList& interpolationCells = e.interpolationCells();
784  const labelListList& cellStencil = e.cellStencil();
785 
787 
788  label nInvalidInterpolated = 0;
789  label nLocal = 0;
790  label nMixed = 0;
791  label nRemote = 0;
792  forAll(interpolationCells, i)
793  {
794  label celli = interpolationCells[i];
795  const labelList& slots = cellStencil[celli];
796 
797  if (slots.empty())
798  {
799  nInvalidInterpolated++;
800  }
801 
802  bool hasLocal = false;
803  bool hasRemote = false;
804 
805  forAll(slots, sloti)
806  {
807  if (slots[sloti] >= cellTypes.size())
808  {
809  hasRemote = true;
810  }
811  else
812  {
813  hasLocal = true;
814  }
815  }
816 
817  if (hasRemote)
818  {
819  if (!hasLocal)
820  {
821  nRemote++;
822  }
823  else
824  {
825  nMixed++;
826  }
827  }
828  else if (hasLocal)
829  {
830  nLocal++;
831  }
832  }
833  reduce(nLocal, sumOp<label>());
834  reduce(nMixed, sumOp<label>());
835  reduce(nRemote, sumOp<label>());
836  reduce(nInvalidInterpolated, sumOp<label>());
837 
838  os << "Overset analysis : nCells : "
839  << returnReduce(cellTypes.size(), sumOp<label>()) << nl
840  << incrIndent
841  << indent << "calculated : " << nCells[cellCellStencil::CALCULATED]
842  << nl
843  << indent << "interpolated : " << nCells[cellCellStencil::INTERPOLATED]
844  << " (interpolated from local:" << nLocal
845  << " mixed local/remote:" << nMixed
846  << " remote:" << nRemote
847  << " special:" << nInvalidInterpolated << ")" << nl
848  << indent << "hole : " << nCells[cellCellStencil::HOLE] << nl;
849 
850  if (nCells[cellCellStencil::SPECIAL] || nCells[cellCellStencil::POROUS])
851  {
852  os << indent << "special : " << nCells[cellCellStencil::SPECIAL]
853  << nl
854  << indent << "porous : " << nCells[cellCellStencil::POROUS]
855  << nl;
856  }
857  os << decrIndent << endl;
858 
859  return os;
860 }
861 
862 
863 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
virtual const wordHashSet & nonInterpolatedFields() const
Return the names of any (stencil or mesh specific) fields that.
dictionary dict
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
word findInstance(const fileName &dir, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null) const
Return time instance (location) of dir that contains the file name (eg, used in reading mesh data)...
Definition: Time.C:725
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
void seedCell(const label cellI, const scalar wantedFraction, bitSet &isFront, scalarField &fraction) const
Seed faces of cell with wantedFraction (if higher than current)
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static autoPtr< cellCellStencil > New(const fvMesh &, const dictionary &dict, const bool update=true)
New function which constructs and returns pointer to a.
void setUpFront(const labelList &allCellTypes, bitSet &isFront) const
Set up front using allCellTypes.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:432
const cellList & cells() const
virtual ~cellCellStencil()
Destructor.
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:830
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool localStencil(const labelUList &) const
Helper: is stencil fully local.
void setUpFrontOnOversetPatch(const labelList &allCellTypes, bitSet &isFront) const
Set up front on overset patches.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
label findMin(const ListType &input, label start=0)
Linear search for the index of the min element, similar to std::min_element but for lists and returns...
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
bool ends_with(char c) const
True if string ends with given character (cf. C++20)
Definition: string.H:460
static labelList count(const label size, const labelUList &lst)
Count occurrences (in parallel)
void suppressMotionFields()
Helper: populate nonInterpolatedFields_ with motion solver.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
A class for handling words, derived from Foam::string.
Definition: word.H:63
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const labelList & cellTypes
Definition: setCellMask.H:27
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static word baseName(const word &name)
Helper: strip off trailing _0.
virtual bool update()=0
Update stencils. Return false if nothing changed.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
mesh update()
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
static void globalCellCells(const globalIndex &gi, const polyMesh &mesh, const boolList &isValidDonor, const labelList &selectedCells, labelListList &cellCells, pointListList &cellCellCentres)
Helper: create cell-cell addressing in global numbering.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
vector point
Point is a vector.
Definition: point.H:37
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
label nCells() const noexcept
Number of mesh cells.
static const Enum< cellType > cellTypeNames_
Mode type names.
Nothing to be read.
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.
messageStream Info
Information stream (stdout output on master, null elsewhere)
static void listCombineGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
List< point > pointList
List of point.
Definition: pointList.H:32
void walkFront(const globalIndex &globalCells, const scalar layerRelax, const labelListList &allStencil, labelList &allCellTypes, scalarField &allWeight, const scalarList &compactCellVol, const labelListList &compactStencil, const labelList &zoneID, const label holeLayers, const label useLayer) const
Surround holes with layer(s) of interpolated cells.
List< label > labelList
A List of labels.
Definition: List.H:62
List< pointList > pointListList
List of pointList.
Definition: pointList.H:35
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
Request registration (bool: true)
List< bool > boolList
A List of bools.
Definition: List.H:60
Do not request registration (bool: false)
Namespace for OpenFOAM.
static tmp< volScalarField > createField(const fvMesh &mesh, const scalar val)
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