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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  subsetMesh
26 
27 Description
28  Selects a section of mesh based on a cellSet.
29 
30  The utility sub-sets the mesh to choose only a part of interest. Check
31  the setSet/cellSet/topoSet utilities to see how to select cells based on
32  various shapes.
33 
34  The mesh will subset all points, faces and cells needed to make a sub-mesh
35  but will not preserve attached boundary types.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "fvMeshSubset.H"
40 #include "argList.H"
41 #include "cellSet.H"
42 #include "IOobjectList.H"
43 #include "volFields.H"
44 #include "topoDistanceData.H"
45 #include "FaceCellWave.H"
46 
47 using namespace Foam;
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 labelList nearestPatch(const polyMesh& mesh, const labelList& patchIDs)
52 {
53  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
54 
55  // Count number of faces in exposedPatchIDs
56  label nFaces = 0;
57  forAll(patchIDs, i)
58  {
59  const polyPatch& pp = pbm[patchIDs[i]];
60  nFaces += pp.size();
61  }
62 
63  // Field on cells and faces.
66 
67  // Start of changes
68  labelList patchFaces(nFaces);
69  List<topoDistanceData> patchData(nFaces);
70  nFaces = 0;
71  forAll(patchIDs, i)
72  {
73  label patchI = patchIDs[i];
74  const polyPatch& pp = pbm[patchI];
75 
76  forAll(pp, i)
77  {
78  patchFaces[nFaces] = pp.start()+i;
79  patchData[nFaces] = topoDistanceData(patchI, 0);
80  nFaces++;
81  }
82  }
83 
84  // Propagate information inwards
86  (
87  mesh,
88  patchFaces,
89  patchData,
90  faceData,
91  cellData,
93  );
94 
95  // And extract
96 
97  labelList nearest(mesh.nFaces());
98 
99  bool haveWarned = false;
100  forAll(faceData, faceI)
101  {
102  if (!faceData[faceI].valid(deltaCalc.data()))
103  {
104  if (!haveWarned)
105  {
107  << "Did not visit some faces, e.g. face " << faceI
108  << " at " << mesh.faceCentres()[faceI] << endl
109  << "Using patch " << patchIDs[0] << " as nearest"
110  << endl;
111  haveWarned = true;
112  }
113  nearest[faceI] = patchIDs[0];
114  }
115  else
116  {
117  nearest[faceI] = faceData[faceI].data();
118  }
119  }
120 
121  return nearest;
122 }
123 
124 
125 template<class Type>
126 void subsetVolFields
127 (
128  const fvMeshSubset& subsetter,
129  const wordList& fieldNames,
131 )
132 {
133  const fvMesh& baseMesh = subsetter.baseMesh();
134 
135  forAll(fieldNames, i)
136  {
137  const word& fieldName = fieldNames[i];
138 
139  Info<< "Subsetting field " << fieldName << endl;
140 
142  (
143  IOobject
144  (
145  fieldName,
146  baseMesh.time().timeName(),
147  baseMesh,
150  ),
151  baseMesh
152  );
153 
154  subFields.set(i, subsetter.interpolate(fld));
155  }
156 }
157 
158 
159 template<class Type>
161 (
162  const fvMeshSubset& subsetter,
163  const wordList& fieldNames,
165 )
166 {
167  const fvMesh& baseMesh = subsetter.baseMesh();
168 
169  forAll(fieldNames, i)
170  {
171  const word& fieldName = fieldNames[i];
172 
173  Info<< "Subsetting field " << fieldName << endl;
174 
176  (
177  IOobject
178  (
179  fieldName,
180  baseMesh.time().timeName(),
181  baseMesh,
184  ),
185  baseMesh
186  );
187 
188  subFields.set(i, subsetter.interpolate(fld));
189  }
190 }
191 
192 
193 template<class Type>
195 (
196  const fvMeshSubset& subsetter,
197  const pointMesh& pMesh,
198  const wordList& fieldNames,
200 )
201 {
202  const fvMesh& baseMesh = subsetter.baseMesh();
203 
204  forAll(fieldNames, i)
205  {
206  const word& fieldName = fieldNames[i];
207 
208  Info<< "Subsetting field " << fieldName << endl;
209 
211  (
212  IOobject
213  (
214  fieldName,
215  baseMesh.time().timeName(),
216  baseMesh,
219  ),
220  pMesh
221  );
222 
223  subFields.set(i, subsetter.interpolate(fld));
224  }
225 }
226 
227 
228 template<class Type>
230 (
231  const fvMeshSubset& subsetter,
232  const wordList& fieldNames,
234 )
235 {
236  const fvMesh& baseMesh = subsetter.baseMesh();
237 
238  forAll(fieldNames, i)
239  {
240  const word& fieldName = fieldNames[i];
241 
242  Info<< "Subsetting field " << fieldName << endl;
243 
245  (
246  IOobject
247  (
248  fieldName,
249  baseMesh.time().timeName(),
250  baseMesh,
253  ),
254  baseMesh
255  );
256 
257  subFields.set(i, subsetter.interpolate(fld));
258  }
259 }
260 
261 
262 
263 int main(int argc, char *argv[])
264 {
266  (
267  "select a mesh subset based on a cellSet"
268  );
269 
270  #include "addOverwriteOption.H"
271  #include "addRegionOption.H"
272  argList::validArgs.append("cellSet");
274  (
275  "patch",
276  "name",
277  "add exposed internal faces to specified patch instead of to "
278  "'oldInternalFaces'"
279  );
281  (
282  "patches",
283  "names",
284  "add exposed internal faces to nearest of specified patches"
285  " instead of to 'oldInternalFaces'"
286  );
288  (
289  "resultTime",
290  "time",
291  "specify a time for the resulting mesh"
292  );
293  #include "setRootCase.H"
294  #include "createTime.H"
295  runTime.functionObjects().off();
296 
297  Foam::word meshRegionName = polyMesh::defaultRegion;
298  args.optionReadIfPresent("region", meshRegionName);
299 
300  #include "createNamedMesh.H"
301 
302 
303  const word setName = args[1];
304 
305  word meshInstance = mesh.pointsInstance();
306  word fieldsInstance = runTime.timeName();
307 
308  const bool overwrite = args.optionFound("overwrite");
309  const bool specifiedInstance = args.optionReadIfPresent
310  (
311  "resultTime",
312  fieldsInstance
313  );
314  if (specifiedInstance)
315  {
316  // Set both mesh and field to this time
317  meshInstance = fieldsInstance;
318  }
319 
320 
321  Info<< "Reading cell set from " << setName << endl << endl;
322 
323  // Create mesh subsetting engine
324  fvMeshSubset subsetter(mesh);
325 
326  labelList exposedPatchIDs;
327 
328  if (args.optionFound("patch"))
329  {
330  const word patchName = args["patch"];
331 
332  exposedPatchIDs = labelList
333  (
334  1,
335  mesh.boundaryMesh().findPatchID(patchName)
336  );
337 
338  if (exposedPatchIDs[0] == -1)
339  {
341  << nl << "Valid patches are " << mesh.boundaryMesh().names()
342  << exit(FatalError);
343  }
344 
345  Info<< "Adding exposed internal faces to patch " << patchName << endl
346  << endl;
347  }
348  else if (args.optionFound("patches"))
349  {
350  const wordReList patchNames(args.optionRead<wordReList>("patches"));
351 
352  exposedPatchIDs = mesh.boundaryMesh().patchSet(patchNames).sortedToc();
353 
354  Info<< "Adding exposed internal faces to nearest of patches "
355  << patchNames << endl << endl;
356  }
357  else
358  {
359  Info<< "Adding exposed internal faces to a patch called"
360  << " \"oldInternalFaces\" (created if necessary)" << endl
361  << endl;
362  exposedPatchIDs = labelList(1, label(-1));
363  }
364 
365 
366  cellSet currentSet(mesh, setName);
367 
368  if (exposedPatchIDs.size() == 1)
369  {
370  subsetter.setLargeCellSubset(currentSet, exposedPatchIDs[0], true);
371  }
372  else
373  {
374 
375  // Find per face the nearest patch
376  labelList nearestExposedPatch(nearestPatch(mesh, exposedPatchIDs));
377 
378  labelList region(mesh.nCells(), 0);
379  forAllConstIter(cellSet, currentSet, iter)
380  {
381  region[iter.key()] = 1;
382  }
383 
384  labelList exposedFaces(subsetter.getExposedFaces(region, 1, true));
385  subsetter.setLargeCellSubset
386  (
387  region,
388  1,
389  exposedFaces,
390  UIndirectList<label>(nearestExposedPatch, exposedFaces)(),
391  true
392  );
393  }
394 
395 
396  IOobjectList objects(mesh, runTime.timeName());
397 
398 
399  // Read vol fields and subset
400  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
401 
402  wordList scalarNames(objects.names(volScalarField::typeName));
403  PtrList<volScalarField> scalarFlds(scalarNames.size());
404  subsetVolFields(subsetter, scalarNames, scalarFlds);
405 
406  wordList vectorNames(objects.names(volVectorField::typeName));
407  PtrList<volVectorField> vectorFlds(vectorNames.size());
408  subsetVolFields(subsetter, vectorNames, vectorFlds);
409 
410  wordList sphericalTensorNames
411  (
412  objects.names(volSphericalTensorField::typeName)
413  );
414  PtrList<volSphericalTensorField> sphericalTensorFlds
415  (
416  sphericalTensorNames.size()
417  );
418  subsetVolFields(subsetter, sphericalTensorNames, sphericalTensorFlds);
419 
420  wordList symmTensorNames(objects.names(volSymmTensorField::typeName));
421  PtrList<volSymmTensorField> symmTensorFlds(symmTensorNames.size());
422  subsetVolFields(subsetter, symmTensorNames, symmTensorFlds);
423 
424  wordList tensorNames(objects.names(volTensorField::typeName));
425  PtrList<volTensorField> tensorFlds(tensorNames.size());
426  subsetVolFields(subsetter, tensorNames, tensorFlds);
427 
428 
429  // Read surface fields and subset
430  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
431 
432  wordList surfScalarNames(objects.names(surfaceScalarField::typeName));
433  PtrList<surfaceScalarField> surfScalarFlds(surfScalarNames.size());
434  subsetSurfaceFields(subsetter, surfScalarNames, surfScalarFlds);
435 
436  wordList surfVectorNames(objects.names(surfaceVectorField::typeName));
437  PtrList<surfaceVectorField> surfVectorFlds(surfVectorNames.size());
438  subsetSurfaceFields(subsetter, surfVectorNames, surfVectorFlds);
439 
440  wordList surfSphericalTensorNames
441  (
442  objects.names(surfaceSphericalTensorField::typeName)
443  );
444  PtrList<surfaceSphericalTensorField> surfSphericalTensorFlds
445  (
446  surfSphericalTensorNames.size()
447  );
449  (
450  subsetter,
451  surfSphericalTensorNames,
452  surfSphericalTensorFlds
453  );
454 
455  wordList surfSymmTensorNames
456  (
457  objects.names(surfaceSymmTensorField::typeName)
458  );
459  PtrList<surfaceSymmTensorField> surfSymmTensorFlds
460  (
461  surfSymmTensorNames.size()
462  );
463  subsetSurfaceFields(subsetter, surfSymmTensorNames, surfSymmTensorFlds);
464 
465  wordList surfTensorNames(objects.names(surfaceTensorField::typeName));
466  PtrList<surfaceTensorField> surfTensorFlds(surfTensorNames.size());
467  subsetSurfaceFields(subsetter, surfTensorNames, surfTensorFlds);
468 
469 
470  // Read point fields and subset
471  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
472 
473  const pointMesh& pMesh = pointMesh::New(mesh);
474 
475  wordList pointScalarNames(objects.names(pointScalarField::typeName));
476  PtrList<pointScalarField> pointScalarFlds(pointScalarNames.size());
477  subsetPointFields(subsetter, pMesh, pointScalarNames, pointScalarFlds);
478 
479  wordList pointVectorNames(objects.names(pointVectorField::typeName));
480  PtrList<pointVectorField> pointVectorFlds(pointVectorNames.size());
481  subsetPointFields(subsetter, pMesh, pointVectorNames, pointVectorFlds);
482 
483  wordList pointSphericalTensorNames
484  (
485  objects.names(pointSphericalTensorField::typeName)
486  );
487  PtrList<pointSphericalTensorField> pointSphericalTensorFlds
488  (
489  pointSphericalTensorNames.size()
490  );
492  (
493  subsetter,
494  pMesh,
495  pointSphericalTensorNames,
496  pointSphericalTensorFlds
497  );
498 
499  wordList pointSymmTensorNames
500  (
501  objects.names(pointSymmTensorField::typeName)
502  );
503  PtrList<pointSymmTensorField> pointSymmTensorFlds
504  (
505  pointSymmTensorNames.size()
506  );
508  (
509  subsetter,
510  pMesh,
511  pointSymmTensorNames,
512  pointSymmTensorFlds
513  );
514 
515  wordList pointTensorNames(objects.names(pointTensorField::typeName));
516  PtrList<pointTensorField> pointTensorFlds(pointTensorNames.size());
517  subsetPointFields(subsetter, pMesh, pointTensorNames, pointTensorFlds);
518 
519 
520  // Read dimensioned fields and subset
521  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
522 
523  typedef volScalarField::DimensionedInternalField dimScalType;
524  wordList scalarDimNames(objects.names(dimScalType::typeName));
525  PtrList<dimScalType> scalarDimFlds(scalarDimNames.size());
526  subsetDimensionedFields(subsetter, scalarDimNames, scalarDimFlds);
527 
528  typedef volVectorField::DimensionedInternalField dimVecType;
529  wordList vectorDimNames(objects.names(dimVecType::typeName));
530  PtrList<dimVecType> vectorDimFlds(vectorDimNames.size());
531  subsetDimensionedFields(subsetter, vectorDimNames, vectorDimFlds);
532 
534  wordList sphericalTensorDimNames(objects.names(dimSphereType::typeName));
535  PtrList<dimSphereType> sphericalTensorDimFlds
536  (
537  sphericalTensorDimNames.size()
538  );
540  (
541  subsetter,
542  sphericalTensorDimNames,
543  sphericalTensorDimFlds
544  );
545 
546  typedef volSymmTensorField::DimensionedInternalField dimSymmTensorType;
547  wordList symmTensorDimNames(objects.names(dimSymmTensorType::typeName));
548  PtrList<dimSymmTensorType> symmTensorDimFlds(symmTensorDimNames.size());
549  subsetDimensionedFields(subsetter, symmTensorDimNames, symmTensorDimFlds);
550 
551  typedef volTensorField::DimensionedInternalField dimTensorType;
552  wordList tensorDimNames(objects.names(dimTensorType::typeName));
553  PtrList<dimTensorType> tensorDimFlds(tensorDimNames.size());
554  subsetDimensionedFields(subsetter, tensorDimNames, tensorDimFlds);
555 
556 
557  // Write mesh and fields to new time
558  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
559 
560  if (overwrite || specifiedInstance)
561  {
562  runTime.setTime(instant(fieldsInstance), 0);
563  subsetter.subMesh().setInstance(meshInstance);
564  }
565  else
566  {
567  runTime++;
568  }
569 
570  Info<< "Writing subsetted mesh and fields to time " << runTime.timeName()
571  << endl;
572  subsetter.subMesh().write();
573 
574 
575  // Subsetting adds 'subset' prefix. Rename field to be like original.
576  forAll(scalarFlds, i)
577  {
578  scalarFlds[i].rename(scalarNames[i]);
579  scalarFlds[i].write();
580  }
581  forAll(vectorFlds, i)
582  {
583  vectorFlds[i].rename(vectorNames[i]);
584  vectorFlds[i].write();
585  }
586  forAll(sphericalTensorFlds, i)
587  {
588  sphericalTensorFlds[i].rename(sphericalTensorNames[i]);
589  sphericalTensorFlds[i].write();
590  }
591  forAll(symmTensorFlds, i)
592  {
593  symmTensorFlds[i].rename(symmTensorNames[i]);
594  symmTensorFlds[i].write();
595  }
596  forAll(tensorFlds, i)
597  {
598  tensorFlds[i].rename(tensorNames[i]);
599  tensorFlds[i].write();
600  }
601 
602  // Surface ones.
603  forAll(surfScalarFlds, i)
604  {
605  surfScalarFlds[i].rename(surfScalarNames[i]);
606  surfScalarFlds[i].write();
607  }
608  forAll(surfVectorFlds, i)
609  {
610  surfVectorFlds[i].rename(surfVectorNames[i]);
611  surfVectorFlds[i].write();
612  }
613  forAll(surfSphericalTensorFlds, i)
614  {
615  surfSphericalTensorFlds[i].rename(surfSphericalTensorNames[i]);
616  surfSphericalTensorFlds[i].write();
617  }
618  forAll(surfSymmTensorFlds, i)
619  {
620  surfSymmTensorFlds[i].rename(surfSymmTensorNames[i]);
621  surfSymmTensorFlds[i].write();
622  }
623  forAll(surfTensorNames, i)
624  {
625  surfTensorFlds[i].rename(surfTensorNames[i]);
626  surfTensorFlds[i].write();
627  }
628 
629  // Point ones
630  forAll(pointScalarFlds, i)
631  {
632  pointScalarFlds[i].rename(pointScalarNames[i]);
633  pointScalarFlds[i].write();
634  }
635  forAll(pointVectorFlds, i)
636  {
637  pointVectorFlds[i].rename(pointVectorNames[i]);
638  pointVectorFlds[i].write();
639  }
640  forAll(pointSphericalTensorFlds, i)
641  {
642  pointSphericalTensorFlds[i].rename(pointSphericalTensorNames[i]);
643  pointSphericalTensorFlds[i].write();
644  }
645  forAll(pointSymmTensorFlds, i)
646  {
647  pointSymmTensorFlds[i].rename(pointSymmTensorNames[i]);
648  pointSymmTensorFlds[i].write();
649  }
650  forAll(pointTensorNames, i)
651  {
652  pointTensorFlds[i].rename(pointTensorNames[i]);
653  pointTensorFlds[i].write();
654  }
655 
656  // DimensionedFields
657  forAll(scalarDimFlds, i)
658  {
659  scalarDimFlds[i].rename(scalarDimNames[i]);
660  scalarDimFlds[i].write();
661  }
662  forAll(vectorDimFlds, i)
663  {
664  vectorDimFlds[i].rename(vectorDimNames[i]);
665  vectorDimFlds[i].write();
666  }
667  forAll(sphericalTensorDimFlds, i)
668  {
669  sphericalTensorDimFlds[i].rename(sphericalTensorDimNames[i]);
670  sphericalTensorDimFlds[i].write();
671  }
672  forAll(symmTensorDimFlds, i)
673  {
674  symmTensorDimFlds[i].rename(symmTensorDimNames[i]);
675  symmTensorDimFlds[i].write();
676  }
677  forAll(tensorDimFlds, i)
678  {
679  tensorDimFlds[i].rename(tensorDimNames[i]);
680  tensorDimFlds[i].write();
681  }
682 
683 
684  Info<< "\nEnd\n" << endl;
685 
686  return 0;
687 }
688 
689 
690 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
Foam::fvMeshSubset::setLargeCellSubset
void setLargeCellSubset(const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
Set the subset from all cells with region == currentRegion.
Definition: fvMeshSubset.C:795
volFields.H
Foam::fvMeshSubset::interpolate
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Definition: fvMeshSubsetInterpolate.C:43
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
topoDistanceData.H
Foam::topoDistanceData
For use with FaceCellWave. Determines topological distance to starting faces.
Definition: topoDistanceData.H:53
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::argList::addOption
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:108
Foam::argList::addNote
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:139
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
fvMeshSubset.H
addOverwriteOption.H
subsetVolFields
void subsetVolFields(const fvMeshSubset &subsetter, const wordList &fieldNames, PtrList< GeometricField< Type, fvPatchField, volMesh > > &subFields)
Definition: subsetMesh.C:124
Foam::globalMeshData::nTotalCells
label nTotalCells() const
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:394
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:527
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::nl
static const char nl
Definition: Ostream.H:260
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:66
argList.H
addRegionOption.H
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
patchNames
wordList patchNames(nPatches)
Foam::fvMeshSubset::getExposedFaces
labelList getExposedFaces(const labelList &region, const label currentRegion, const bool syncCouples=true) const
Two step subsetting.
Definition: fvMeshSubset.C:1405
Foam::FatalError
error FatalError
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:48
Foam::HashTable::sortedToc
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:216
subsetSurfaceFields
void subsetSurfaceFields(const fvMeshSubset &subsetter, const wordList &fieldNames, PtrList< GeometricField< Type, fvsPatchField, surfaceMesh > > &subFields)
Definition: subsetMesh.C:158
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
main
int main(int argc, char *argv[])
Definition: subsetMesh.C:260
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::FaceCellWave
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:75
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
nearestPatch
labelList nearestPatch(const polyMesh &mesh, const labelList &patchIDs)
Definition: subsetMesh.C:48
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::fieldNames
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
Definition: ReadFields.C:36
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:130
FaceCellWave.H
Foam::argList::optionRead
T optionRead(const word &opt) const
Read a value from the named option.
Definition: argListI.H:187
Foam::fvMeshSubset::subMesh
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubset.C:1472
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:64
cellSet.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::fvMeshSubset::baseMesh
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:216
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:750
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
Foam::argList::optionReadIfPresent
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
subsetDimensionedFields
void subsetDimensionedFields(const fvMeshSubset &subsetter, const wordList &fieldNames, PtrList< DimensionedField< Type, volMesh > > &subFields)
Definition: subsetMesh.C:227
subsetPointFields
void subsetPointFields(const fvMeshSubset &subsetter, const pointMesh &pMesh, const wordList &fieldNames, PtrList< GeometricField< Type, pointPatchField, pointMesh > > &subFields)
Definition: subsetMesh.C:192
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51