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-2024 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 "fvMeshSubsetter.H" // Not fvMeshSubset (need two-step subsetting)
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 {
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 DimensionedField/GeometricField
153 //
154 template<class FieldType>
155 PtrList<FieldType> subsetFields
156 (
157  const fvMeshSubset& subsetter,
158  const IOobjectList& objects
159 )
160 {
161  const fvMesh& baseMesh = subsetter.baseMesh();
162 
163  const UPtrList<const IOobject> fieldObjects
164  (
165  objects.csorted<FieldType>()
166  );
167 
168  PtrList<FieldType> subFields(fieldObjects.size());
169 
170  label nFields = 0;
171  for (const IOobject& io : fieldObjects)
172  {
173  if (!nFields)
174  {
175  Info<< "Subsetting " << FieldType::typeName << " (";
176  }
177  else
178  {
179  Info<< ' ';
180  }
181  Info<< io.name();
182 
183  FieldType fld
184  (
185  IOobject
186  (
187  io.name(),
188  baseMesh.time().timeName(),
189  baseMesh,
193  ),
194  baseMesh
195  );
196 
197  subFields.set(nFields, subsetter.interpolate(fld));
198  auto& subField = subFields[nFields];
199  ++nFields;
200 
201  // Subsetting adds 'subset' prefix - rename to match original.
202  subField.rename(io.name());
203  }
204 
205  if (nFields)
206  {
207  Info<< ')' << nl;
208  }
209 
210  return subFields;
211 }
212 
213 
214 // Subset point fields
215 template<class FieldType>
216 PtrList<FieldType> subsetFields
217 (
218  const fvMeshSubset& subsetter,
219  const IOobjectList& objects,
220  const pointMesh& pMesh
221 )
222 {
223  //const fvMesh& baseMesh = subsetter.baseMesh();
224 
225  const UPtrList<const IOobject> fieldObjects
226  (
227  objects.csorted<FieldType>()
228  );
229 
230  PtrList<FieldType> subFields(fieldObjects.size());
231 
232  label nFields = 0;
233  for (const IOobject& io : fieldObjects)
234  {
235  if (!nFields)
236  {
237  Info<< "Subsetting " << FieldType::typeName << " (";
238  }
239  else
240  {
241  Info<< ' ';
242  }
243  Info<< io.name();
244 
245  FieldType fld
246  (
247  IOobject
248  (
249  io.name(),
250  pMesh.thisDb().time().timeName(),
251  pMesh.thisDb(),
255  ),
256  pMesh
257  );
258 
259  subFields.set(nFields, subsetter.interpolate(fld));
260  auto& subField = subFields[nFields];
261  ++nFields;
262 
263  // Subsetting adds 'subset' prefix - rename to match original.
264  subField.rename(io.name());
265  }
266 
267  if (nFields)
268  {
269  Info<< ')' << nl;
270  }
271 
272  return subFields;
273 }
274 
275 
276 template<class TopoSet>
277 void subsetTopoSets
278 (
279  const fvMesh& mesh,
280  const IOobjectList& objects,
281  const labelList& map,
282  const fvMesh& subMesh,
283  PtrList<TopoSet>& subSets
284 )
285 {
286  // Read original sets
287  PtrList<TopoSet> sets;
288  ReadFields<TopoSet>(objects, sets);
289 
290  subSets.resize_null(sets.size());
291 
292  forAll(sets, seti)
293  {
294  const TopoSet& set = sets[seti];
295 
296  Info<< "Subsetting " << set.type() << " " << set.name() << endl;
297 
299  subset.reserve(Foam::min(set.size(), map.size()));
300 
301  // Map the data
302  forAll(map, i)
303  {
304  if (set.found(map[i]))
305  {
306  subset.insert(i);
307  }
308  }
309 
310  subSets.set
311  (
312  seti,
313  new TopoSet
314  (
315  subMesh,
316  set.name(),
317  std::move(subset),
319  )
320  );
321  }
322 }
323 
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326 
327 int main(int argc, char *argv[])
328 {
330  (
331  "Create a mesh subset for a particular region of interest based on a"
332  " cellSet or cellZone(s) specified as the first command argument.\n"
333  "See setSet/topoSet utilities on how to select cells based on"
334  " various shapes."
335  );
336 
337  #include "addOverwriteOption.H"
338  #include "addRegionOption.H"
340  (
341  "cell-selection",
342  "The cellSet name, but with the -zone option this is interpreted"
343  " to be a cellZone selection by name(s) or regex.\n"
344  "Eg 'mixer' or '( mixer \"moving.*\" )'"
345  );
346 
348  (
349  "patch",
350  "name",
351  "Add exposed internal faces to specified patch"
352  " instead of \"oldInternalFaces\""
353  );
355  (
356  "patches",
357  "wordRes",
358  "Add exposed internal faces to closest of specified patches"
359  " instead of \"oldInternalFaces\""
360  );
362  (
363  "exclude-patches",
364  "wordRes",
365  "Exclude single or multiple patches from the -patches selection"
366  );
368  (
369  "zone",
370  "Subset with cellZone(s) instead of cellSet."
371  " The command argument may be a list of words or regexs"
372  );
374  (
375  "resultTime",
376  "time",
377  "Specify a time for the resulting mesh"
378  );
379 
380  argList::noFunctionObjects(); // Never use function objects
381 
382  #include "setRootCase.H"
383  #include "createTime.H"
384 
385  #include "createNamedMesh.H"
386  // Make sure pointMesh gets constructed/read as well
388 
389  // arg[1] = word (cellSet) or wordRes (cellZone)
390  // const word selectionName = args[1];
391 
392  word meshInstance = mesh.pointsInstance();
393  word fieldsInstance = runTime.timeName();
394 
395  const bool useCellZone = args.found("zone");
396  const bool overwrite = args.found("overwrite");
397  const bool specifiedInstance = args.readIfPresent
398  (
399  "resultTime",
400  fieldsInstance
401  );
402  if (specifiedInstance)
403  {
404  // Set both mesh and field to this time
405  meshInstance = fieldsInstance;
406  }
407 
408 
409  // Default exposed patch id
410  labelList exposedPatchIDs(one{}, -1);
411 
412  wordRes includePatches, excludePatches;
413 
414  if (!args.readListIfPresent<wordRe>("patches", includePatches))
415  {
416  if (args.found("patch"))
417  {
418  includePatches.resize(1);
419  includePatches.front() = args.get<word>("patch");
420  }
421  }
422  args.readListIfPresent<wordRe>("exclude-patches", excludePatches);
423 
424  if (includePatches.size() == 1 && includePatches.front().isLiteral())
425  {
426  // Select a single patch - no exclude possible
427  exposedPatchIDs.front() =
428  getExposedPatchId(mesh, includePatches.front());
429  }
430  else if (!includePatches.empty())
431  {
432  // Patches selected (sorted order)
433  exposedPatchIDs =
434  mesh.boundaryMesh().indices(includePatches, excludePatches);
435 
436  // Only retain initial, non-processor patches
437  const label nNonProcessor
438  (
440  );
441 
442  forAll(exposedPatchIDs, i)
443  {
444  if (exposedPatchIDs[i] > nNonProcessor)
445  {
446  exposedPatchIDs.resize(i);
447  break;
448  }
449  }
450 
451  const wordList allPatchNames(mesh.boundaryMesh().names());
452 
453  Info<< "Adding exposed internal faces to nearest of patches:" << nl
454  << " include: " << flatOutput(includePatches) << nl
455  << " exclude: " << flatOutput(excludePatches) << nl
456  << nl;
457 
458  if (exposedPatchIDs.empty())
459  {
461  << nl << "No patches matched. Patches: "
462  << flatOutput(allPatchNames) << nl
463  << exit(FatalError);
464  }
465  }
466  else
467  {
468  Info<< "Adding exposed internal faces to patch \""
470  << "\" (created if necessary)" << nl
471  << nl;
472  }
473 
474 
475  autoPtr<cellSet> cellSetPtr;
476 
477  // arg[1] can be a word (for cellSet) or wordRes (for cellZone)
478 
479  wordRes zoneNames;
480  if (useCellZone)
481  {
482  wordRes selectionNames(args.getList<wordRe>(1));
483  zoneNames.transfer(selectionNames);
484 
485  Info<< "Using cellZone " << flatOutput(zoneNames) << nl << endl;
486 
487  if (mesh.cellZones().findIndex(zoneNames) == -1)
488  {
490  << "No cellZones found: " << flatOutput(zoneNames) << nl << nl
491  << exit(FatalError);
492  }
493  }
494  else
495  {
496  const word selectionName = args[1];
497 
498  Info<< "Using cellSet " << selectionName << nl << endl;
499 
500  cellSetPtr.emplace(mesh, selectionName);
501  }
502 
503 
504  // Two-step mesh subsetting engine
505  fvMeshSubsetter subsetter(mesh);
506 
507  {
508  bitSet selectedCells =
509  (
510  cellSetPtr
511  ? BitSetOps::create(mesh.nCells(), *cellSetPtr)
512  : mesh.cellZones().selection(zoneNames)
513  );
514 
515  if (exposedPatchIDs.size() == 1)
516  {
517  // Single patch for exposed faces (syncPar)
518  subsetter.reset(selectedCells, exposedPatchIDs.front(), true);
519  }
520  else
521  {
522  // The nearest patch per face
523  labelList nearestExposedPatch(nearestPatch(mesh, exposedPatchIDs));
524 
525  labelList exposedFaces
526  (
527  subsetter.getExposedFaces(selectedCells, true) // syncPar
528  );
529 
530  subsetter.setCellSubset
531  (
532  selectedCells,
533  exposedFaces,
534  labelUIndList(nearestExposedPatch, exposedFaces)(),
535  true // syncPar
536  );
537  }
538 
539  FixedList<label, 2> cellCount;
540  cellCount[0] = subsetter.subMesh().nCells();
541  cellCount[1] = mesh.nCells();
542  reduce(cellCount, sumOp<label>());
543 
544  Info<< "Subset " << cellCount[0] << " of " << cellCount[1]
545  << " cells" << nl << nl;
546  }
547 
548 
549  IOobjectList objects(mesh, runTime.timeName());
550 
551 
552  // Read fields and subset
553  #undef createSubsetFields
554  #define createSubsetFields(FieldType, Variable) \
555  PtrList<FieldType> Variable \
556  ( \
557  subsetFields<FieldType>(subsetter, objects) \
558  );
559 
560 
561  // Read vol fields and subset
562  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
563  createSubsetFields(volScalarField, vScalarFlds);
564  createSubsetFields(volVectorField, vVectorFlds);
565  createSubsetFields(volSphericalTensorField, vSphTensorFlds);
566  createSubsetFields(volSymmTensorField, vSymmTensorFlds);
567  createSubsetFields(volTensorField, vTensorFlds);
568 
569  // Read surface fields and subset
570  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
571  createSubsetFields(surfaceScalarField, sScalarFlds);
572  createSubsetFields(surfaceVectorField, sVectorFlds);
573  createSubsetFields(surfaceSphericalTensorField, sSphTensorFlds);
574  createSubsetFields(surfaceSymmTensorField, sSymmTensorFlds);
575  createSubsetFields(surfaceTensorField, sTensorFlds);
576 
577  // Read dimensioned fields and subset
578  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
579  createSubsetFields(volScalarField::Internal, dScalarFlds);
580  createSubsetFields(volVectorField::Internal, dVectorFlds);
581  createSubsetFields(volSphericalTensorField::Internal, dSphTensorFlds);
582  createSubsetFields(volSymmTensorField::Internal, dSymmTensorFlds);
583  createSubsetFields(volTensorField::Internal, dTensorFlds);
584 
585 
586  // Read point fields and subset
587  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
588 
590 
591  #undef createSubsetFields
592  #define createSubsetFields(FieldType, Variable) \
593  PtrList<FieldType> Variable \
594  ( \
595  subsetFields<FieldType>(subsetter, objects, pMesh) \
596  );
597 
598  createSubsetFields(pointScalarField, pScalarFlds);
599  createSubsetFields(pointVectorField, pVectorFlds);
600  createSubsetFields(pointSphericalTensorField, pSphTensorFlds);
601  createSubsetFields(pointSymmTensorField, pSymmTensorFlds);
602  createSubsetFields(pointTensorField, pTensorFlds);
603 
604  #undef createSubsetFields
605 
606 
607  // Read topoSets and subset
608  // ~~~~~~~~~~~~~~~~~~~~~~~~
609 
610  PtrList<cellSet> cellSets;
611  PtrList<faceSet> faceSets;
612  PtrList<pointSet> pointSets;
613 
614  {
615  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
616  if (cellSetPtr)
617  {
618  objects.remove(*cellSetPtr);
619  }
620  subsetTopoSets
621  (
622  mesh,
623  objects,
624  subsetter.cellMap(),
625  subsetter.subMesh(),
626  cellSets
627  );
628  subsetTopoSets
629  (
630  mesh,
631  objects,
632  subsetter.faceMap(),
633  subsetter.subMesh(),
634  faceSets
635  );
636  subsetTopoSets
637  (
638  mesh,
639  objects,
640  subsetter.pointMap(),
641  subsetter.subMesh(),
642  pointSets
643  );
644  }
645 
646 
647  // Write mesh and fields to new time
648  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
649 
650  if (overwrite || specifiedInstance)
651  {
652  runTime.setTime(instant(fieldsInstance), 0);
653  subsetter.subMesh().setInstance(meshInstance);
654  topoSet::setInstance(meshInstance, cellSets);
655  topoSet::setInstance(meshInstance, faceSets);
656  topoSet::setInstance(meshInstance, pointSets);
657  }
658  else
659  {
660  ++runTime;
661  subsetter.subMesh().setInstance(runTime.timeName());
662  }
663 
664  Info<< "Writing subsetted mesh and fields to time " << runTime.timeName()
665  << endl;
666  subsetter.subMesh().write();
668 
669  auto* subPointMeshPtr =
670  subsetter.subMesh().thisDb().findObject<pointMesh>
671  (
672  pointMesh::typeName
673  );
674  if (subPointMeshPtr)
675  {
676  pointMesh& subPointMesh = const_cast<pointMesh&>(*subPointMeshPtr);
677  subPointMesh.setInstance(subsetter.subMesh().facesInstance());
678  subPointMesh.write();
679  }
680 
681 
682  // Volume fields
683  for (const auto& fld : vScalarFlds) { fld.write(); }
684  for (const auto& fld : vVectorFlds) { fld.write(); }
685  for (const auto& fld : vSphTensorFlds) { fld.write(); }
686  for (const auto& fld : vSymmTensorFlds) { fld.write(); }
687  for (const auto& fld : vTensorFlds) { fld.write(); }
688 
689  // Surface fields
690  for (const auto& fld : sScalarFlds) { fld.write(); }
691  for (const auto& fld : sVectorFlds) { fld.write(); }
692  for (const auto& fld : sSphTensorFlds) { fld.write(); }
693  for (const auto& fld : sSymmTensorFlds) { fld.write(); }
694  for (const auto& fld : sTensorFlds) { fld.write(); }
695 
696  // Dimensioned fields
697  for (const auto& fld : dScalarFlds) { fld.write(); }
698  for (const auto& fld : dVectorFlds) { fld.write(); }
699  for (const auto& fld : dSphTensorFlds) { fld.write(); }
700  for (const auto& fld : dSymmTensorFlds) { fld.write(); }
701  for (const auto& fld : dTensorFlds) { fld.write(); }
702 
703  // Point fields
704  for (const auto& fld : pScalarFlds) { fld.write(); }
705  for (const auto& fld : pVectorFlds) { fld.write(); }
706  for (const auto& fld : pSphTensorFlds) { fld.write(); }
707  for (const auto& fld : pSymmTensorFlds) { fld.write(); }
708  for (const auto& fld : pTensorFlds) { fld.write(); }
709 
710  Info<< "\nEnd\n" << endl;
711 
712  return 0;
713 }
714 
715 
716 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
label patchId(-1)
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const polyBoundaryMesh & pbm
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:561
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:476
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
static word exposedPatchName
Name for exposed internal faces (default: oldInternalFaces)
Definition: fvMeshSubset.H:173
T * front()
The first entry in the list.
Definition: UILList.H:141
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
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:844
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:150
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:337
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:600
label findIndex(const wordRe &key) const
Zone index for the first match, return -1 if not found.
Definition: ZoneMesh.C:695
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
void setCellSubset(const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset()
Definition: fvMeshSubset.H:389
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:265
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:206
List< T > getList(const label index) const
Get a List of values from the argument at index.
Required Classes.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:388
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Field reading functions for post-processing utilities.
const fvMesh & baseMesh() const noexcept
Original mesh.
Definition: fvMeshSubsetI.H:23
T & emplace(Args &&... args)
Reset with emplace construction. Return reference to the new content.
Definition: autoPtrI.H:59
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
label nFaces() const noexcept
Number of mesh faces.
const labelList & faceMap() const
Return face map.
Definition: fvMeshSubsetI.H:65
void reduce(T &value, [[maybe_unused]] BinaryOp bop, [[maybe_unused]] const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce)
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
labelList indices(const wordRe &matcher, const bool useGroups=true) const
The (sorted) patch indices for all matches, optionally matching patch groups.
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubsetI.H:41
void resize_null(const label newLen)
Set the addressed list to the given size, deleting all existing entries. Afterwards the list contains...
Definition: PtrListI.H:96
Required Classes.
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:838
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
wordList names() const
Return a list of patch names.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:201
Reading is optional [identical to LAZY_READ].
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
const labelList & pointMap() const
Return point map.
Definition: fvMeshSubsetI.H:57
void reset()
Reset subMesh and all maps. Same as clear()
Definition: fvMeshSubset.C:491
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:399
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1296
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
const labelList & cellMap() const
Return cell map.
Definition: fvMeshSubsetI.H:84
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings...
Definition: wordRe.H:78
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1068
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
bool readListIfPresent(const word &optName, List< T > &list) const
If named option is present, get a List of values treating a single entry like a list of size 1...
Definition: argListI.H:387
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:75
label nTotalCells() const noexcept
Total global number of mesh cells.
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))
const vectorField & faceCentres() const
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name...
Definition: instant.H:53
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: pointMesh.C:129
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
label nCells() const noexcept
Number of mesh cells.
autoPtr< IOobject > remove(const IOobject &io)
Remove object from the list by its IOobject::name().
Definition: IOobjectList.H:302
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Automatically write from objectRegistry::writeObject()
For use with FaceCellWave. Determines topological distance to starting faces. Templated on passive tr...
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:365
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
messageStream Info
Information stream (stdout output on master, null elsewhere)
bitSet selection(const labelUList &zoneIds) const
Return all elements (cells, faces, points) contained in the listed zones.
Definition: ZoneMesh.C:800
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
bitSet create(const label n, const labelHashSet &locations, const bool on=true)
Create a bitSet with length n with the specified on locations.
Definition: BitOps.C:228
label nNonProcessor() const
The number of patches before the first processor patch.
Foam::argList args(argc, argv)
Extends Foam::fvMeshSubset with two-step subsetting (uses polyTopoChange modification).
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
bool found
Do not request registration (bool: false)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:56
static void setInstance(const fileName &instance, Container &items)
Helper: set instance on all items in container.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225