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-2022 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  if (!mesh.foundObject<labelIOList>("zoneID"))
116  {
117  labelIOList* zoneIDPtr = new labelIOList
118  (
119  IOobject
120  (
121  "zoneID",
124  mesh,
127  ),
128  mesh.nCells()
129  );
130  labelIOList& zoneID = *zoneIDPtr;
131 
132  volScalarField volZoneID
133  (
134  IOobject
135  (
136  "zoneID",
137  mesh.time().findInstance(mesh.dbDir(), "zoneID"),
138  mesh,
141  false
142  ),
143  mesh
144  );
145  forAll(volZoneID, celli)
146  {
147  zoneID[celli] = label(volZoneID[celli]);
148  }
149 
150  zoneIDPtr->store();
151  }
152  return mesh.lookupObject<labelIOList>("zoneID");
153 }
154 
155 
157 (
158  const label size,
159  const labelUList& lst
160 )
161 {
162  labelList count(size, Zero);
163  forAll(lst, i)
164  {
165  count[lst[i]]++;
166  }
168  return count;
169 }
170 
173 {
174  return nonInterpolatedFields_;
175 }
176 
179 {
180  return nonInterpolatedFields_;
181 }
182 
183 
184 bool Foam::cellCellStencil::localStencil(const labelUList& slots) const
185 {
186  forAll(slots, i)
187  {
188  if (slots[i] >= mesh_.nCells())
189  {
190  return false;
191  }
192  }
193  return true;
194 }
195 
196 
198 (
199  const globalIndex& gi,
200  const polyMesh& mesh,
201  const boolList& isValidCell,
202  const labelList& selectedCells,
203  labelListList& cellCells,
204  pointListList& cellCellCentres
205 )
206 {
207  // For selected cells determine the face neighbours (in global numbering)
208 
209  const pointField& cellCentres = mesh.cellCentres();
210  const labelList& faceOwner = mesh.faceOwner();
211  const labelList& faceNeighbour = mesh.faceNeighbour();
212  const cellList& cells = mesh.cells();
213 
214 
215  // 1. Determine global cell number on other side of coupled patches
216 
217  labelList globalCellIDs(identity(gi.localSize(), gi.localStart()));
218 
219  labelList nbrGlobalCellIDs;
221  (
222  mesh,
223  globalCellIDs,
224  nbrGlobalCellIDs
225  );
226  pointField nbrCellCentres;
228  (
229  mesh,
230  cellCentres,
231  nbrCellCentres
232  );
233 
234  boolList nbrIsValidCell;
236  (
237  mesh,
238  isValidCell,
239  nbrIsValidCell
240  );
241 
242 
243  // 2. Collect cell and all its neighbours
244 
245  cellCells.setSize(mesh.nCells());
246  cellCellCentres.setSize(cellCells.size());
247 
248  forAll(selectedCells, i)
249  {
250  label celli = selectedCells[i];
251 
252  const cell& cFaces = cells[celli];
253  labelList& stencil = cellCells[celli];
254  pointList& stencilPoints = cellCellCentres[celli];
255  stencil.setSize(cFaces.size()+1);
256  stencilPoints.setSize(stencil.size());
257  label compacti = 0;
258 
259  // First entry is cell itself
260  if (isValidCell[celli])
261  {
262  stencil[compacti] = globalCellIDs[celli];
263  stencilPoints[compacti++] = cellCentres[celli];
264  }
265 
266  // Other entries are cell neighbours
267  forAll(cFaces, i)
268  {
269  label facei = cFaces[i];
270  label bFacei = facei-mesh.nInternalFaces();
271  label own = faceOwner[facei];
272  label nbrCelli;
273  point nbrCc;
274  bool isValid = false;
275  if (bFacei >= 0)
276  {
277  nbrCelli = nbrGlobalCellIDs[bFacei];
278  nbrCc = nbrCellCentres[bFacei];
279  isValid = nbrIsValidCell[bFacei];
280  }
281  else
282  {
283  if (own != celli)
284  {
285  nbrCelli = gi.toGlobal(own);
286  nbrCc = cellCentres[own];
287  isValid = isValidCell[own];
288  }
289  else
290  {
291  label nei = faceNeighbour[facei];
292  nbrCelli = gi.toGlobal(nei);
293  nbrCc = cellCentres[nei];
294  isValid = isValidCell[nei];
295  }
296  }
297 
298  if (isValid)
299  {
300  SubList<label> current(stencil, compacti);
301  if (!current.found(nbrCelli))
302  {
303  stencil[compacti] = nbrCelli;
304  stencilPoints[compacti++] = nbrCc;
305  }
306  }
307  }
308  stencil.setSize(compacti);
309  stencilPoints.setSize(compacti);
310  }
311 }
312 
314 (
315  const label celli,
316  const scalar wantedFraction,
317  bitSet& isFront,
318  scalarField& fraction
319 ) const
320 {
321  const cell& cFaces = mesh_.cells()[celli];
322  forAll(cFaces, i)
323  {
324  const label nbrFacei = cFaces[i];
325  if (fraction[nbrFacei] < wantedFraction)
326  {
327  fraction[nbrFacei] = wantedFraction;
328  isFront.set(nbrFacei);
329  }
330  }
331 }
332 
333 
335 (
336  const labelList& allCellTypes,
337  bitSet& isFront
338 ) const
339 {
340  const labelList& own = mesh_.faceOwner();
341  const labelList& nei = mesh_.faceNeighbour();
342 
343  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
344  {
345  const label ownType = allCellTypes[own[facei]];
346  const label neiType = allCellTypes[nei[facei]];
347  if
348  (
349  (ownType == HOLE && neiType != HOLE)
350  || (ownType != HOLE && neiType == HOLE)
351  )
352  {
353  isFront.set(facei);
354  }
355  }
356 
357  labelList nbrCellTypes;
358  syncTools::swapBoundaryCellList(mesh_, allCellTypes, nbrCellTypes);
359 
360  for
361  (
362  label facei = mesh_.nInternalFaces();
363  facei < mesh_.nFaces();
364  facei++
365  )
366  {
367  const label ownType = allCellTypes[own[facei]];
368  const label neiType = nbrCellTypes[facei-mesh_.nInternalFaces()];
369 
370  if
371  (
372  (ownType == HOLE && neiType != HOLE)
373  || (ownType != HOLE && neiType == HOLE)
374  )
375  {
376  isFront.set(facei);
377  }
378  }
379 }
380 
381 
383 (
384  const labelList& allCellTypes,
385  bitSet& isFront
386 ) const
387 {
388  const fvBoundaryMesh& fvm = mesh_.boundary();
389  // 'overset' patches
390  forAll(fvm, patchi)
391  {
392  if (isA<oversetFvPatch>(fvm[patchi]))
393  {
394  const labelList& fc = fvm[patchi].faceCells();
395  forAll(fc, i)
396  {
397  const label celli = fc[i];
398  if (allCellTypes[celli] == INTERPOLATED)
399  {
400  // Note that acceptors might have been marked hole if
401  // there are no donors in which case we do not want to
402  // walk this out. This is an extreme situation.
403  isFront.set(fvm[patchi].start()+i);
404  }
405  }
406  }
407  }
408 }
409 
410 
412 (
413  const globalIndex& globalCells,
414  const label layerRelax,
415  const labelListList& allStencil,
416  labelList& allCellTypes,
417  scalarField& allWeight,
418  const scalarList& compactCellVol,
419  const labelListList& compactStencil,
420  const labelList& zoneID,
421  const label holeLayers,
422  const label useLayer
423 ) const
424 {
425  if (useLayer > holeLayers)
426  {
427  FatalErrorInFunction<< "useLayer: " << useLayer
428  << " is larger than : " << holeLayers
429  << abort(FatalError);
430  }
431 
432  if (useLayer == 0)
433  {
434  FatalErrorInFunction<< "useLayer: " << useLayer
435  << " can not be zero."
436  << abort(FatalError);
437  }
438 
439  // Current front
440  bitSet isFront(mesh_.nFaces());
441 
442  // List of cellTYpes
443  DynamicList<labelList> dataSet;
444  // List of cellTypes for increasing layers
445  DynamicList<scalarField> dAlllWeight;
446  // List of average volumen ration interpolated/donor
447  DynamicList<scalar> daverageVolRatio;
448 
449  // Counting holes of different layers
450  DynamicList<label> dHoles;
451 
452  const scalarField& V = mesh_.V();
453 
454  // Current layer
455  labelField nLayer(allCellTypes.size(), 0);
456 
457  for (label currLayer = 1; currLayer < holeLayers+1; ++currLayer)
458  {
459  // Set up original front
460  isFront = false;
461 
462  setUpFront(allCellTypes, isFront);
463 
464  labelList allCellTypesWork(allCellTypes);
465 
466  bitSet isFrontWork(isFront);
467  label nCurrLayer = currLayer;
468 
469  while (nCurrLayer > 1 && returnReduce(isFrontWork.any(), orOp<bool>()))
470  {
471  bitSet newIsFront(mesh_.nFaces());
472  forAll(isFrontWork, facei)
473  {
474  if (isFrontWork.test(facei))
475  {
476  const label own = mesh_.faceOwner()[facei];
477 
478  if (allCellTypesWork[own] != HOLE)
479  {
480  allCellTypesWork[own] = HOLE;
481  newIsFront.set(mesh_.cells()[own]);
482  }
483 
484  if (mesh_.isInternalFace(facei))
485  {
486  const label nei = mesh_.faceNeighbour()[facei];
487  if (allCellTypesWork[nei] != HOLE)
488  {
489  allCellTypesWork[nei] = HOLE;
490  newIsFront.set(mesh_.cells()[nei]);
491  }
492  }
493  }
494  }
495  syncTools::syncFaceList(mesh_, newIsFront, orEqOp<unsigned int>());
496 
497  isFrontWork.transfer(newIsFront);
498 
499  nCurrLayer--;
500  }
501 
502 
503  if ((debug&2) && (mesh_.time().outputTime()))
504  {
505  tmp<volScalarField> tfld
506  (
507  createField
508  (
509  mesh_,
510  "allCellTypesWork_Holes" + name(currLayer),
511  allCellTypesWork
512  )
513  );
514  tfld().write();
515  }
516 
517  if (currLayer == 1)
518  {
519  setUpFrontOnOversetPatch(allCellTypes, isFront);
520  }
521 
522  if (currLayer > 1)
523  {
524  isFront = false;
525  setUpFrontOnOversetPatch(allCellTypesWork, isFront);
526  setUpFront(allCellTypesWork, isFront);
527  }
528 
529  // Current interpolation fraction
530  scalarField fraction(mesh_.nFaces(), Zero);
531 
532  // Ratio between Inter/donor
533  scalarField volRatio(allCellTypes.size(), Zero);
534 
535  forAll(isFront, facei)
536  {
537  if (isFront.test(facei))
538  {
539  fraction[facei] = 1.0;
540  }
541  }
542 
543  scalarField allWeightWork(allCellTypes.size(), Zero);
544  bitSet nHoles(allCellTypes.size());
545 
546  while (returnReduce(isFront.any(), orOp<bool>()))
547  {
548  // Interpolate cells on front
549  bitSet newIsFront(mesh_.nFaces());
550  scalarField newFraction(fraction);
551 
552  forAll(isFront, facei)
553  {
554  if (isFront.test(facei))
555  {
556  const label own = mesh_.faceOwner()[facei];
557 
558  if (allCellTypesWork[own] != HOLE)
559  {
560  if (allWeightWork[own] < fraction[facei])
561  {
562  // Cell wants to become interpolated (if sufficient
563  // stencil, otherwise becomes hole)
564  if (allStencil[own].size())
565  {
566  nLayer[own] = currLayer;
567 
568  allWeightWork[own] = fraction[facei];
569  allCellTypesWork[own] = INTERPOLATED;
570 
571  const label donorId = compactStencil[own][0];
572 
573  volRatio[own] = V[own]/compactCellVol[donorId];
574 
575  seedCell
576  (
577  own,
578  fraction[facei]-layerRelax,
579  newIsFront,
580  newFraction
581  );
582 
583  nHoles[own] = true;
584  }
585  else
586  {
587  allWeightWork[own] = 0.0;
588  allCellTypesWork[own] = HOLE;
589  // Add faces of cell as new front
590  seedCell
591  (
592  own,
593  1.0,
594  newIsFront,
595  newFraction
596  );
597  }
598  }
599  }
600  if (mesh_.isInternalFace(facei))
601  {
602  label nei = mesh_.faceNeighbour()[facei];
603  if (allCellTypesWork[nei] != HOLE)
604  {
605  if (allWeightWork[nei] < fraction[facei])
606  {
607  if (allStencil[nei].size())
608  {
609  nLayer[nei] = currLayer;
610 
611  allWeightWork[nei] = fraction[facei];
612  allCellTypesWork[nei] = INTERPOLATED;
613 
614  const label donorId = compactStencil[nei][0];
615 
616  volRatio[nei] = V[nei]/compactCellVol[donorId];
617 
618  seedCell
619  (
620  nei,
621  fraction[facei]-layerRelax,
622  newIsFront,
623  newFraction
624  );
625  }
626  else
627  {
628  allWeightWork[nei] = 0.0;
629  allCellTypesWork[nei] = HOLE;
630  nHoles[nei] = true;
631  seedCell
632  (
633  nei,
634  1.0,
635  newIsFront,
636  newFraction
637  );
638  }
639  }
640  }
641  }
642  }
643  }
644 
645  syncTools::syncFaceList(mesh_, newIsFront, orEqOp<unsigned int>());
646  syncTools::syncFaceList(mesh_, newFraction, maxEqOp<scalar>());
647 
648  isFront.transfer(newIsFront);
649  fraction.transfer(newFraction);
650  }
651 
652  if ((debug&2) && (mesh_.time().outputTime()))
653  {
654  tmp<volScalarField> tfld
655  (
656  createField
657  (
658  mesh_,
659  "allCellTypesWork_Layers" + name(currLayer),
660  allCellTypesWork
661  )
662  );
663  tfld().write();
664  }
665 
666  dAlllWeight.append(allWeightWork);
667  dataSet.append(allCellTypesWork);
668 
669  // Counting interpolated cells
670  label count = 1;
671  forAll (volRatio, i)
672  {
673  if (volRatio[i] > 0)
674  {
675  count++;
676  }
677  }
678  label nCount(count);
679  reduce(nCount, sumOp<label>());
680 
681  const scalar aveVol = mag(gSum(volRatio));
682 
683  daverageVolRatio.append(aveVol/nCount);
684 
685  // Check holes number. A sudden increase occurs when the walk leaks
686  // out from the obstacle
687  label nTotalHoles(nHoles.count());
688  reduce(nTotalHoles, sumOp<label>());
689  dHoles.append(nTotalHoles);
690 
691  // Check the increase between this layer and the last one
692  // if over 50% breaks the layer loop.
693  if
694  (
695  (currLayer > 1)
696  & (nTotalHoles > 2.0*dHoles[currLayer - 1])
697  )
698  {
699  break;
700  }
701  }
702 
703 
704  if ((debug&2) && (mesh_.time().outputTime()))
705  {
706  tmp<volScalarField> tfld
707  (
708  createField(mesh_, "walkFront_layers", nLayer)
709  );
710  tfld().write();
711  }
712 
713  // Try to find the best averageVolRatio the further from the initial set
714  // As this one is next to HOLES
715  scalarList averageVolRatio;
716  averageVolRatio.transfer(daverageVolRatio);
717  label minVolId = findMin(averageVolRatio);
718 
719  if (useLayer > -1)
720  {
721  Info<< nl << " Number of layers : " << averageVolRatio.size() << nl
722  << " Average volumetric ratio : " << averageVolRatio << nl
723  << " Number of holes cells : " << dHoles << nl
724  << " Using layer : " << useLayer << nl << endl;
725  }
726 
727  if (useLayer != -1)
728  {
729  minVolId = useLayer - 1;
730  }
731 
732  allCellTypes.transfer(dataSet[minVolId]);
733  allWeight.transfer(dAlllWeight[minVolId]);
734 }
735 
736 
737 // * * * * * * * * * * * * * Ostream operator * * * * * * * * * * * * * * * //
738 
739 template<>
740 Foam::Ostream& Foam::operator<<
741 (
742  Ostream& os,
743  const InfoProxy<cellCellStencil>& ic
744 )
745 {
746  const cellCellStencil& e = ic.t_;
747 
748  const labelUList& cellTypes = e.cellTypes();
749  const labelUList& interpolationCells = e.interpolationCells();
750  const labelListList& cellStencil = e.cellStencil();
751 
753 
754  label nInvalidInterpolated = 0;
755  label nLocal = 0;
756  label nMixed = 0;
757  label nRemote = 0;
758  forAll(interpolationCells, i)
759  {
760  label celli = interpolationCells[i];
761  const labelList& slots = cellStencil[celli];
762 
763  if (slots.empty())
764  {
765  nInvalidInterpolated++;
766  }
767 
768  bool hasLocal = false;
769  bool hasRemote = false;
770 
771  forAll(slots, sloti)
772  {
773  if (slots[sloti] >= cellTypes.size())
774  {
775  hasRemote = true;
776  }
777  else
778  {
779  hasLocal = true;
780  }
781  }
782 
783  if (hasRemote)
784  {
785  if (!hasLocal)
786  {
787  nRemote++;
788  }
789  else
790  {
791  nMixed++;
792  }
793  }
794  else if (hasLocal)
795  {
796  nLocal++;
797  }
798  }
799  reduce(nLocal, sumOp<label>());
800  reduce(nMixed, sumOp<label>());
801  reduce(nRemote, sumOp<label>());
802  reduce(nInvalidInterpolated, sumOp<label>());
803 
804  Info<< "Overset analysis : nCells : "
805  << returnReduce(cellTypes.size(), sumOp<label>()) << nl
806  << incrIndent
807  << indent << "calculated : " << nCells[cellCellStencil::CALCULATED]
808  << nl
809  << indent << "interpolated : " << nCells[cellCellStencil::INTERPOLATED]
810  << " (interpolated from local:" << nLocal
811  << " mixed local/remote:" << nMixed
812  << " remote:" << nRemote
813  << " special:" << nInvalidInterpolated << ")" << nl
814  << indent << "hole : " << nCells[cellCellStencil::HOLE] << nl
815  << indent << "special : " << nCells[cellCellStencil::SPECIAL] << nl
816  << decrIndent << endl;
817 
818  return os;
819 }
820 
821 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
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:118
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
void walkFront(const globalIndex &globalCells, const label 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.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:449
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:853
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:777
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:439
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:578
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1110
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:463
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
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:402
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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:824
Ignore writing from objectRegistry::writeObject()
List< point > pointList
A List of points.
Definition: pointList.H:61
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
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.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:80
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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:84
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:439
static labelList count(const label size, const labelUList &lst)
Count occurrences (in parallel)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
const cellShapeList & cells
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:61
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1104
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:55
static word baseName(const word &name)
Helper: strip off trailing _0.
virtual bool update()=0
Update stencils. Return false if nothing changed.
List< pointList > pointListList
A List of pointList.
Definition: pointList.H:63
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:467
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...
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)
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
static void listCombineGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
List< label > labelList
A List of labels.
Definition: List.H:62
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:458
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:615
List< bool > boolList
A List of bools.
Definition: List.H:60
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:41
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:38
Namespace for OpenFOAM.
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:157