subsetMesh.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  subsetMesh
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Create a mesh subset for a particular region of interest based on a
35  cellSet or cellZone.
36 
37  See setSet/topoSet utilities on how to define select cells based on
38  various shapes.
39 
40  Will subset all points, faces and cells needed to make a sub-mesh,
41  but not preserve attached boundary types.
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #include "fvMeshSubset.H"
46 #include "argList.H"
47 #include "IOobjectList.H"
48 #include "volFields.H"
49 #include "topoDistanceData.H"
50 #include "FaceCellWave.H"
51 #include "cellSet.H"
52 #include "faceSet.H"
53 #include "pointSet.H"
54 #include "ReadFields.H"
55 #include "processorMeshes.H"
56 
57 using namespace Foam;
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 // Get the exposed patchId or define the exposedPatchName in fvMeshSubset
62 label getExposedPatchId(const polyMesh& mesh, const word& patchName)
63 {
64  const label patchId = mesh.boundaryMesh().findPatchID(patchName);
65 
66  if (patchId == -1)
67  {
69  }
70 
71  Info<< "Adding exposed internal faces to "
72  << (patchId == -1 ? "new" : "existing")
73  << " patch \"" << patchName << "\"" << nl << endl;
74 
75  return patchId;
76 }
77 
78 
79 labelList nearestPatch(const polyMesh& mesh, const labelList& patchIDs)
80 {
81  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
82 
83  // Count number of faces in exposedPatchIDs
84  label nFaces = 0;
85  for (const label patchi : patchIDs)
86  {
87  nFaces += pbm[patchi].size();
88  }
89 
90  // Field on cells and faces.
93 
94  // Start of changes
95  labelList patchFaces(nFaces);
96  List<topoDistanceData<label>> patchData(nFaces);
97  nFaces = 0;
98  for (const label patchi : patchIDs)
99  {
100  const polyPatch& pp = pbm[patchi];
101 
102  forAll(pp, i)
103  {
104  patchFaces[nFaces] = pp.start()+i;
105  patchData[nFaces] = topoDistanceData<label>(0, patchi);
106  ++nFaces;
107  }
108  }
109 
110  // Propagate information inwards
112  (
113  mesh,
114  patchFaces,
115  patchData,
116  faceData,
117  cellData,
119  );
120 
121  // And extract
122 
123  labelList nearest(mesh.nFaces());
124 
125  bool haveWarned = false;
126  forAll(faceData, faceI)
127  {
128  if (!faceData[faceI].valid(deltaCalc.data()))
129  {
130  if (!haveWarned)
131  {
133  << "Did not visit some faces, e.g. face " << faceI
134  << " at " << mesh.faceCentres()[faceI] << nl
135  << "Using patch " << patchIDs[0] << " as nearest"
136  << endl;
137  haveWarned = true;
138  }
139  nearest[faceI] = patchIDs[0];
140  }
141  else
142  {
143  nearest[faceI] = faceData[faceI].data();
144  }
145  }
146 
147  return nearest;
148 }
149 
150 
151 //
152 // Subset field-type, availability information cached
153 // in the availableFields hashtable.
154 //
155 template<class Type, template<class> class PatchField, class GeoMesh>
156 void subsetFields
157 (
158  const fvMeshSubset& subsetter,
159  HashTable<wordHashSet>& availableFields,
161 )
162 {
164  const word fieldType = FieldType::typeName;
165 
166  const wordList fieldNames = availableFields(fieldType).sortedToc();
167  subFields.setSize(fieldNames.size());
168 
169  const fvMesh& baseMesh = subsetter.baseMesh();
170 
171  label nFields = 0;
172  for (const word& fieldName : fieldNames)
173  {
174  if (!nFields)
175  {
176  Info<< "Subsetting " << fieldType << " (";
177  }
178  else
179  {
180  Info<< ' ';
181  }
182  Info<< fieldName;
183 
184  FieldType fld
185  (
186  IOobject
187  (
188  fieldName,
189  baseMesh.time().timeName(),
190  baseMesh,
193  ),
194  baseMesh
195  );
196 
197  subFields.set(nFields, subsetter.interpolate(fld));
198 
199  // Subsetting adds 'subset' prefix - rename to match original.
200  subFields[nFields].rename(fieldName);
201 
202  ++nFields;
203  }
204 
205  if (nFields)
206  {
207  Info<< ')' << nl;
208  }
209 }
210 
211 
212 template<class Type>
213 void subsetPointFields
214 (
215  const fvMeshSubset& subsetter,
216  const pointMesh& pMesh,
217  HashTable<wordHashSet>& availableFields,
219 )
220 {
222  const word fieldType = FieldType::typeName;
223 
224  const wordList fieldNames = availableFields(fieldType).sortedToc();
225  subFields.setSize(fieldNames.size());
226 
227  const fvMesh& baseMesh = subsetter.baseMesh();
228 
229  label nFields = 0;
230  for (const word& fieldName : fieldNames)
231  {
232  if (!nFields)
233  {
234  Info<< "Subsetting " << fieldType << " (";
235  }
236  else
237  {
238  Info<< ' ';
239  }
240  Info<< fieldName;
241 
242  FieldType fld
243  (
244  IOobject
245  (
246  fieldName,
247  baseMesh.time().timeName(),
248  baseMesh,
251  ),
252  pMesh
253  );
254 
255  subFields.set(nFields, subsetter.interpolate(fld));
256 
257  // Subsetting adds 'subset' prefix - rename to match original.
258  subFields[nFields].rename(fieldName);
259 
260  ++nFields;
261  }
262 
263  if (nFields)
264  {
265  Info<< ')' << nl;
266  }
267 }
268 
269 
270 template<class Type>
271 void subsetDimensionedFields
272 (
273  const fvMeshSubset& subsetter,
274  HashTable<wordHashSet>& availableFields,
276 )
277 {
278  typedef DimensionedField<Type, volMesh> FieldType;
279  const word fieldType = FieldType::typeName;
280 
281  const wordList fieldNames = availableFields(fieldType).sortedToc();
282  subFields.setSize(fieldNames.size());
283 
284  const fvMesh& baseMesh = subsetter.baseMesh();
285 
286  label nFields = 0;
287  for (const word& fieldName : fieldNames)
288  {
289  if (!nFields)
290  {
291  Info<< "Subsetting " << fieldType << " (";
292  }
293  else
294  {
295  Info<< ' ';
296  }
297  Info<< fieldName;
298 
299  FieldType fld
300  (
301  IOobject
302  (
303  fieldName,
304  baseMesh.time().timeName(),
305  baseMesh,
308  ),
309  baseMesh
310  );
311 
312  subFields.set(nFields, subsetter.interpolate(fld));
313 
314  // Subsetting adds 'subset' prefix - rename to match original.
315  subFields[nFields].rename(fieldName);
316 
317  ++nFields;
318  }
319 
320  if (nFields)
321  {
322  Info<< ')' << nl;
323  }
324 }
325 
326 
327 template<class TopoSet>
328 void subsetTopoSets
329 (
330  const fvMesh& mesh,
331  const IOobjectList& objects,
332  const labelList& map,
333  const fvMesh& subMesh,
334  PtrList<TopoSet>& subSets
335 )
336 {
337  // Read original sets
338  PtrList<TopoSet> sets;
339  ReadFields<TopoSet>(objects, sets);
340 
341  subSets.setSize(sets.size());
342  forAll(sets, seti)
343  {
344  const TopoSet& set = sets[seti];
345 
346  Info<< "Subsetting " << set.type() << " " << set.name() << endl;
347 
348  labelHashSet subset(2*min(set.size(), map.size()));
349 
350  forAll(map, i)
351  {
352  if (set.found(map[i]))
353  {
354  subset.insert(i);
355  }
356  }
357 
358  subSets.set
359  (
360  seti,
361  new TopoSet
362  (
363  subMesh,
364  set.name(),
365  std::move(subset),
367  )
368  );
369  }
370 }
371 
372 
373 
374 int main(int argc, char *argv[])
375 {
377  (
378  "Create a mesh subset for a particular region of interest based on a"
379  " cellSet or cellZone(s) specified as the first command argument.\n"
380  "See setSet/topoSet utilities on how to select cells based on"
381  " various shapes."
382  );
383 
384  #include "addOverwriteOption.H"
385  #include "addRegionOption.H"
387  (
388  "cell-selection",
389  "The cellSet name, but with the -zone option this is interpreted"
390  " to be a cellZone selection by name(s) or regex.\n"
391  "Eg 'mixer' or '( mixer \"moving.*\" )'"
392  );
393 
395  (
396  "patch",
397  "name",
398  "Add exposed internal faces to specified patch"
399  " instead of \"oldInternalFaces\""
400  );
402  (
403  "patches",
404  "wordRes",
405  "Add exposed internal faces to closest of specified patches"
406  " instead of \"oldInternalFaces\""
407  );
409  (
410  "zone",
411  "Subset with cellZone(s) instead of cellSet."
412  " The command argument may be a list of words or regexs"
413  );
415  (
416  "resultTime",
417  "time",
418  "Specify a time for the resulting mesh"
419  );
420 
421  argList::noFunctionObjects(); // Never use function objects
422 
423  #include "setRootCase.H"
424  #include "createTime.H"
425 
426  #include "createNamedMesh.H"
427 
428  // arg[1] = word (cellSet) or wordRes (cellZone)
429  // const word selectionName = args[1];
430 
431  word meshInstance = mesh.pointsInstance();
432  word fieldsInstance = runTime.timeName();
433 
434  const bool useCellZone = args.found("zone");
435  const bool overwrite = args.found("overwrite");
436  const bool specifiedInstance = args.readIfPresent
437  (
438  "resultTime",
439  fieldsInstance
440  );
441  if (specifiedInstance)
442  {
443  // Set both mesh and field to this time
444  meshInstance = fieldsInstance;
445  }
446 
447 
448  // Default exposed patch id
449  labelList exposedPatchIDs(one{}, -1);
450 
451  if (args.found("patches"))
452  {
453  const wordRes patchNames(args.getList<wordRe>("patches"));
454 
455  if (patchNames.size() == 1 && patchNames.first().isLiteral())
456  {
457  exposedPatchIDs.first() =
458  getExposedPatchId(mesh, patchNames.first());
459  }
460  else
461  {
462  exposedPatchIDs =
463  mesh.boundaryMesh().patchSet(patchNames).sortedToc();
464 
465  Info<< "Adding exposed internal faces to nearest of patches "
466  << flatOutput(patchNames) << nl << endl;
467 
468  if (exposedPatchIDs.empty())
469  {
471  << nl << "No patches matched. Patches: "
472  << mesh.boundaryMesh().names() << nl
473  << exit(FatalError);
474  }
475  }
476  }
477  else if (args.found("patch"))
478  {
479  exposedPatchIDs.first() =
480  getExposedPatchId(mesh, args.get<word>("patch"));
481  }
482  else
483  {
484  Info<< "Adding exposed internal faces to patch \""
486  << "\" (created if necessary)" << nl
487  << nl;
488  }
489 
490 
491  autoPtr<cellSet> cellSetPtr;
492 
493  // arg[1] can be a word (for cellSet) or wordRes (for cellZone)
494 
495  wordRes zoneNames;
496  if (useCellZone)
497  {
498  wordRes selectionNames(args.getList<wordRe>(1));
499  zoneNames.transfer(selectionNames);
500 
501  Info<< "Using cellZone " << flatOutput(zoneNames) << nl << endl;
502 
503  if (mesh.cellZones().findIndex(zoneNames) == -1)
504  {
506  << "No cellZones found: " << flatOutput(zoneNames) << nl << nl
507  << exit(FatalError);
508  }
509  }
510  else
511  {
512  const word selectionName = args[1];
513 
514  Info<< "Using cellSet " << selectionName << nl << endl;
515 
516  cellSetPtr = autoPtr<cellSet>::New(mesh, selectionName);
517  }
518 
519 
520  // Mesh subsetting engine
521  fvMeshSubset subsetter(mesh);
522 
523  {
524  bitSet selectedCells =
525  (
526  cellSetPtr
527  ? BitSetOps::create(mesh.nCells(), *cellSetPtr)
528  : mesh.cellZones().selection(zoneNames)
529  );
530 
531  if (exposedPatchIDs.size() == 1)
532  {
533  // Single patch for exposed faces
534  subsetter.setCellSubset
535  (
536  selectedCells,
537  exposedPatchIDs.first(),
538  true
539  );
540  }
541  else
542  {
543  // The nearest patch per face
544  labelList nearestExposedPatch(nearestPatch(mesh, exposedPatchIDs));
545 
546  labelList exposedFaces
547  (
548  subsetter.getExposedFaces(selectedCells, true)
549  );
550 
551  subsetter.setCellSubset
552  (
553  selectedCells,
554  exposedFaces,
555  labelUIndList(nearestExposedPatch, exposedFaces)(),
556  true
557  );
558  }
559 
560  Info<< "Subset "
561  << returnReduce(subsetter.subMesh().nCells(), sumOp<label>())
562  << " of "
564  << " cells" << nl << nl;
565  }
566 
567 
568  IOobjectList objects(mesh, runTime.timeName());
569  HashTable<wordHashSet> availableFields = objects.classes();
570 
571 
572  // Read vol fields and subset
573  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
574 
575  PtrList<volScalarField> vScalarFlds;
576  subsetFields(subsetter, availableFields, vScalarFlds);
577 
578  PtrList<volVectorField> vVectorFlds;
579  subsetFields(subsetter, availableFields, vVectorFlds);
580 
581  PtrList<volSphericalTensorField> vSphTensorFlds;
582  subsetFields(subsetter, availableFields, vSphTensorFlds);
583 
584  PtrList<volSymmTensorField> vSymmTensorFlds;
585  subsetFields(subsetter, availableFields, vSymmTensorFlds);
586 
587  PtrList<volTensorField> vTensorFlds;
588  subsetFields(subsetter, availableFields, vTensorFlds);
589 
590 
591  // Read surface fields and subset
592  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
593 
594  PtrList<surfaceScalarField> sScalarFlds;
595  subsetFields(subsetter, availableFields, sScalarFlds);
596 
597  PtrList<surfaceVectorField> sVectorFlds;
598  subsetFields(subsetter, availableFields, sVectorFlds);
599 
601  subsetFields(subsetter, availableFields, sSphTensorFlds);
602 
603  PtrList<surfaceSymmTensorField> sSymmTensorFlds;
604  subsetFields(subsetter, availableFields, sSymmTensorFlds);
605 
606  PtrList<surfaceTensorField> sTensorFlds;
607  subsetFields(subsetter, availableFields, sTensorFlds);
608 
609 
610  // Read point fields and subset
611  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
612 
613  const pointMesh& pMesh = pointMesh::New(mesh);
614 
615  PtrList<pointScalarField> pScalarFlds;
616  subsetPointFields(subsetter, pMesh, availableFields, pScalarFlds);
617 
618  PtrList<pointVectorField> pVectorFlds;
619  subsetPointFields(subsetter, pMesh, availableFields, pVectorFlds);
620 
621  PtrList<pointSphericalTensorField> pSphTensorFlds;
622  subsetPointFields(subsetter, pMesh, availableFields, pSphTensorFlds);
623 
624  PtrList<pointSymmTensorField> pSymmTensorFlds;
625  subsetPointFields(subsetter, pMesh, availableFields, pSymmTensorFlds);
626 
627  PtrList<pointTensorField> pTensorFlds;
628  subsetPointFields(subsetter, pMesh, availableFields, pTensorFlds);
629 
630 
631  // Read dimensioned fields and subset
632  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
633 
635  subsetDimensionedFields(subsetter, availableFields, dScalarFlds);
636 
638  subsetDimensionedFields(subsetter, availableFields, dVectorFlds);
639 
641  subsetDimensionedFields(subsetter, availableFields, dSphTensorFlds);
642 
644  subsetDimensionedFields(subsetter, availableFields, dSymmTensorFlds);
645 
647  subsetDimensionedFields(subsetter, availableFields, dTensorFlds);
648 
649 
650  // Read topoSets and subset
651  // ~~~~~~~~~~~~~~~~~~~~~~~~
652 
653  PtrList<cellSet> cellSets;
654  PtrList<faceSet> faceSets;
655  PtrList<pointSet> pointSets;
656 
657  {
658  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
659  if (cellSetPtr)
660  {
661  objects.remove(*cellSetPtr);
662  }
663  subsetTopoSets
664  (
665  mesh,
666  objects,
667  subsetter.cellMap(),
668  subsetter.subMesh(),
669  cellSets
670  );
671  subsetTopoSets
672  (
673  mesh,
674  objects,
675  subsetter.faceMap(),
676  subsetter.subMesh(),
677  faceSets
678  );
679  subsetTopoSets
680  (
681  mesh,
682  objects,
683  subsetter.pointMap(),
684  subsetter.subMesh(),
685  pointSets
686  );
687  }
688 
689 
690  // Write mesh and fields to new time
691  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
692 
693  if (overwrite || specifiedInstance)
694  {
695  runTime.setTime(instant(fieldsInstance), 0);
696  subsetter.subMesh().setInstance(meshInstance);
697  topoSet::setInstance(meshInstance, cellSets);
698  topoSet::setInstance(meshInstance, faceSets);
699  topoSet::setInstance(meshInstance, pointSets);
700  }
701  else
702  {
703  ++runTime;
704  subsetter.subMesh().setInstance(runTime.timeName());
705  }
706 
707  Info<< "Writing subsetted mesh and fields to time " << runTime.timeName()
708  << endl;
709  subsetter.subMesh().write();
711 
712 
713  // Volume fields
714  for (const auto& fld : vScalarFlds) { fld.write(); }
715  for (const auto& fld : vVectorFlds) { fld.write(); }
716  for (const auto& fld : vSphTensorFlds) { fld.write(); }
717  for (const auto& fld : vSymmTensorFlds) { fld.write(); }
718  for (const auto& fld : vTensorFlds) { fld.write(); }
719 
720  // Surface fields
721  for (const auto& fld : sScalarFlds) { fld.write(); }
722  for (const auto& fld : sVectorFlds) { fld.write(); }
723  for (const auto& fld : sSphTensorFlds) { fld.write(); }
724  for (const auto& fld : sSymmTensorFlds) { fld.write(); }
725  for (const auto& fld : sTensorFlds) { fld.write(); }
726 
727  // Point fields
728  for (const auto& fld : pScalarFlds) { fld.write(); }
729  for (const auto& fld : pVectorFlds) { fld.write(); }
730  for (const auto& fld : pSphTensorFlds) { fld.write(); }
731  for (const auto& fld : pSymmTensorFlds) { fld.write(); }
732  for (const auto& fld : pTensorFlds) { fld.write(); }
733 
734  // Dimensioned fields
735  for (const auto& fld : dScalarFlds) { fld.write(); }
736  for (const auto& fld : dVectorFlds) { fld.write(); }
737  for (const auto& fld : dSphTensorFlds) { fld.write(); }
738  for (const auto& fld : dSymmTensorFlds) { fld.write(); }
739  for (const auto& fld : dTensorFlds) { fld.write(); }
740 
741  Info<< "\nEnd\n" << endl;
742 
743  return 0;
744 }
745 
746 
747 // ************************************************************************* //
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::BitSetOps::create
bitSet create(const label n, const labelHashSet &locations, const bool on=true)
Definition: BitOps.C:116
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
topoDistanceData.H
Foam::topoSet::setInstance
static void setInstance(const fileName &instance, Container &)
Definition: topoSetTemplates.C:26
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
Foam::topoDistanceData
For use with FaceCellWave. Determines topological distance to starting faces. Templated on passive tr...
Definition: topoDistanceData.H:49
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Definition: BitOps.C:30
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::ZoneMesh::findIndex
label findIndex(const wordRe &key) const
Definition: ZoneMesh.C:484
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Definition: fvMesh.C:1034
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Foam::fvMeshSubset::cellMap
const labelList & cellMap() const
Definition: fvMeshSubsetI.H:84
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:59
Foam::fvMeshSubset::pointMap
const labelList & pointMap() const
Definition: fvMeshSubsetI.H:57
Foam::fvMeshSubset
Given the original mesh and the list of selected cells, it creates the mesh consisting only of the de...
Definition: fvMeshSubset.H:69
fvMeshSubset.H
addOverwriteOption.H
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Definition: MeshObject.C:41
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::subset
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Foam::fvMeshSubset::exposedPatchName
static word exposedPatchName
Definition: fvMeshSubset.H:157
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Definition: polyMesh.C:845
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:57
Foam::fvMeshSubset::setCellSubset
void setCellSubset(const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
Definition: fvMeshSubset.C:551
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::PtrList::set
const T * set(const label i) const
Definition: PtrList.H:134
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:73
processorMeshes.H
Foam::argList::get
T get(const label index) const
Definition: argListI.H:271
Foam::argList::readIfPresent
bool readIfPresent(const word &optName, T &val) const
Definition: argListI.H:316
Foam::wordRe
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:78
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
Foam::polyBoundaryMesh::names
wordList names() const
Definition: polyBoundaryMesh.C:548
Foam::sumOp
Definition: ops.H:207
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::argList::addArgument
static void addArgument(const string &argName, const string &usage="")
Definition: argList.C:294
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Definition: polyMesh.C:839
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Definition: argList.C:466
Foam::primitiveMesh::nCells
label nCells() const noexcept
Definition: primitiveMeshI.H:89
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:64
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Definition: polyMesh.H:488
argList.H
Foam::IOobjectList::classes
HashTable< wordHashSet > classes() const
Definition: IOobjectList.C:290
Foam::IOobjectList::remove
bool remove(const IOobject &io)
Definition: IOobjectList.C:229
faceSet.H
Foam::List::transfer
void transfer(List< T > &list)
addRegionOption.H
Foam::PtrList::setSize
void setSize(const label newLen)
Definition: PtrList.H:147
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
patchNames
wordList patchNames(nPatches)
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Definition: polyBoundaryMesh.C:758
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::FatalError
error FatalError
createNamedMesh.H
Required Variables.
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam
Definition: atmBoundaryLayer.C:26
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Definition: FlatOutput.H:217
Foam::HashTable::sortedToc
List< Key > sortedToc() const
Definition: HashTable.C:129
Foam::fvMeshSubset::baseMesh
const fvMesh & baseMesh() const noexcept
Definition: fvMeshSubsetI.H:23
Foam::polyPatch::start
label start() const
Definition: polyPatch.H:357
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:101
Foam::argList::addBoolOption
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Definition: argList.C:317
Foam::FaceCellWave
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:74
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::fvMeshSubset::interpolate
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap, const bool allowUnmapped=false)
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
ReadFields.H
Field reading functions for post-processing utilities.
Foam::fvMeshSubset::getExposedFaces
labelList getExposedFaces(const bitSet &selectedCells, const bool syncCouples=true) const
Definition: fvMeshSubset.C:1238
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::Time::setTime
virtual void setTime(const Time &t)
Definition: Time.C:996
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:71
FaceCellWave.H
Foam::processorMeshes::removeFiles
static void removeFiles(const polyMesh &mesh)
Definition: processorMeshes.C:267
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Definition: polyBoundaryMesh.C:857
createTime.H
Foam::argList::getList
List< T > getList(const label index) const
Foam::fvMeshSubset::subMesh
const fvMesh & subMesh() const
Definition: fvMeshSubsetI.H:41
Foam::fvMeshSubset::faceMap
const labelList & faceMap() const
Definition: fvMeshSubsetI.H:65
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::globalMeshData::nTotalCells
label nTotalCells() const noexcept
Definition: globalMeshData.H:367
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Definition: primitiveMeshI.H:83
patchId
label patchId(-1)
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:48
cellSet.H
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Definition: polyMesh.C:1288
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
args
Foam::argList args(argc, argv)
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Definition: polyMeshIO.C:29
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:50
Foam::argList::found
bool found(const word &optName) const
Definition: argListI.H:171
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:54
pointSet.H