subsetMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 the
31  cellSet/topoSet utilities to see how to select cells based on various
32  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 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 
50 template<class Type>
51 void subsetVolFields
52 (
53  const fvMeshSubset& subsetter,
54  const wordList& fieldNames,
55  PtrList<VolField<Type>>& subFields
56 )
57 {
58  const fvMesh& baseMesh = subsetter.baseMesh();
59 
60  forAll(fieldNames, i)
61  {
62  const word& fieldName = fieldNames[i];
63 
64  Info<< "Subsetting field " << fieldName << endl;
65 
67  (
68  IOobject
69  (
70  fieldName,
71  baseMesh.time().name(),
72  baseMesh,
75  ),
76  baseMesh
77  );
78 
79  subFields.set(i, subsetter.interpolate(fld));
80  }
81 }
82 
83 
84 template<class Type>
85 void subsetSurfaceFields
86 (
87  const fvMeshSubset& subsetter,
88  const wordList& fieldNames,
89  PtrList<SurfaceField<Type>>& subFields
90 )
91 {
92  const fvMesh& baseMesh = subsetter.baseMesh();
93 
94  forAll(fieldNames, i)
95  {
96  const word& fieldName = fieldNames[i];
97 
98  Info<< "Subsetting field " << fieldName << endl;
99 
101  (
102  IOobject
103  (
104  fieldName,
105  baseMesh.time().name(),
106  baseMesh,
109  ),
110  baseMesh
111  );
112 
113  subFields.set(i, subsetter.interpolate(fld));
114  }
115 }
116 
117 
118 template<class Type>
119 void subsetPointFields
120 (
121  const fvMeshSubset& subsetter,
122  const pointMesh& pMesh,
123  const wordList& fieldNames,
124  PtrList<PointField<Type>>& subFields
125 )
126 {
127  const fvMesh& baseMesh = subsetter.baseMesh();
128 
129  forAll(fieldNames, i)
130  {
131  const word& fieldName = fieldNames[i];
132 
133  Info<< "Subsetting field " << fieldName << endl;
134 
136  (
137  IOobject
138  (
139  fieldName,
140  baseMesh.time().name(),
141  baseMesh,
144  ),
145  pMesh
146  );
147 
148  subFields.set(i, subsetter.interpolate(fld));
149  }
150 }
151 
152 
153 template<class Type>
154 void subsetDimensionedFields
155 (
156  const fvMeshSubset& subsetter,
157  const wordList& fieldNames,
159 )
160 {
161  const fvMesh& baseMesh = subsetter.baseMesh();
162 
163  forAll(fieldNames, i)
164  {
165  const word& fieldName = fieldNames[i];
166 
167  Info<< "Subsetting field " << fieldName << endl;
168 
170  (
171  IOobject
172  (
173  fieldName,
174  baseMesh.time().name(),
175  baseMesh,
178  ),
179  baseMesh
180  );
181 
182  subFields.set(i, subsetter.interpolate(fld));
183  }
184 }
185 
186 
187 
188 int main(int argc, char *argv[])
189 {
191  (
192  "select a mesh subset based on a cellSet"
193  );
194 
195  #include "addOverwriteOption.H"
196  #include "addRegionOption.H"
197  argList::validArgs.append("cellSet");
199  (
200  "patch",
201  "name",
202  "add exposed internal faces to specified patch instead of to "
203  "'oldInternalFaces'"
204  );
206  (
207  "resultTime",
208  "time",
209  "specify a time for the resulting mesh"
210  );
212  (
213  "noFields",
214  "do not update fields"
215  );
216 
217  #include "setRootCase.H"
219 
220  Foam::word meshRegionName = polyMesh::defaultRegion;
221  args.optionReadIfPresent("region", meshRegionName);
222 
223  #include "createNamedMesh.H"
224 
225 
226  const word setName = args[1];
227 
228  word meshInstance = mesh.pointsInstance();
229  word fieldsInstance = runTime.name();
230 
231  const bool overwrite = args.optionFound("overwrite");
232  const bool specifiedInstance = args.optionReadIfPresent
233  (
234  "resultTime",
235  fieldsInstance
236  );
237  if (specifiedInstance)
238  {
239  // Set both mesh and field to this time
240  meshInstance = fieldsInstance;
241  }
242  const bool fields = !args.optionFound("noFields");
243 
244 
245  Info<< "Reading cell set from " << setName << endl << endl;
246 
247  // Create mesh subsetting engine
248  fvMeshSubset subsetter(mesh);
249 
250  label patchi = -1;
251 
252  if (args.optionFound("patch"))
253  {
254  const word patchName = args["patch"];
255 
256  patchi = mesh.boundaryMesh().findPatchID(patchName);
257 
258  if (patchi == -1)
259  {
261  << nl << "Valid patches are " << mesh.boundaryMesh().names()
262  << exit(FatalError);
263  }
264 
265  Info<< "Adding exposed internal faces to patch " << patchName << endl
266  << endl;
267  }
268  else
269  {
270  Info<< "Adding exposed internal faces to a patch called"
271  << " \"oldInternalFaces\" (created if necessary)" << endl
272  << endl;
273  }
274 
275 
276  cellSet currentSet(mesh, setName);
277  subsetter.setLargeCellSubset(currentSet, patchi, true);
278 
279 
280  if (fields)
281  {
282  IOobjectList objects(mesh, runTime.name());
283 
284  // Read vol fields and subset
285  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
286 
287  wordList scalarNames(objects.names(volScalarField::typeName));
288  PtrList<volScalarField> scalarFlds(scalarNames.size());
289  subsetVolFields(subsetter, scalarNames, scalarFlds);
290 
291  wordList vectorNames(objects.names(volVectorField::typeName));
292  PtrList<volVectorField> vectorFlds(vectorNames.size());
293  subsetVolFields(subsetter, vectorNames, vectorFlds);
294 
295  wordList sphericalTensorNames
296  (
298  );
299  PtrList<volSphericalTensorField> sphericalTensorFlds
300  (
301  sphericalTensorNames.size()
302  );
303  subsetVolFields(subsetter, sphericalTensorNames, sphericalTensorFlds);
304 
305  wordList symmTensorNames(objects.names(volSymmTensorField::typeName));
306  PtrList<volSymmTensorField> symmTensorFlds(symmTensorNames.size());
307  subsetVolFields(subsetter, symmTensorNames, symmTensorFlds);
308 
309  wordList tensorNames(objects.names(volTensorField::typeName));
310  PtrList<volTensorField> tensorFlds(tensorNames.size());
311  subsetVolFields(subsetter, tensorNames, tensorFlds);
312 
313 
314  // Read surface fields and subset
315  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
316 
317  wordList surfScalarNames(objects.names(surfaceScalarField::typeName));
318  PtrList<surfaceScalarField> surfScalarFlds(surfScalarNames.size());
319  subsetSurfaceFields(subsetter, surfScalarNames, surfScalarFlds);
320 
321  wordList surfVectorNames(objects.names(surfaceVectorField::typeName));
322  PtrList<surfaceVectorField> surfVectorFlds(surfVectorNames.size());
323  subsetSurfaceFields(subsetter, surfVectorNames, surfVectorFlds);
324 
325  wordList surfSphericalTensorNames
326  (
328  );
329  PtrList<surfaceSphericalTensorField> surfSphericalTensorFlds
330  (
331  surfSphericalTensorNames.size()
332  );
333  subsetSurfaceFields
334  (
335  subsetter,
336  surfSphericalTensorNames,
337  surfSphericalTensorFlds
338  );
339 
340  wordList surfSymmTensorNames
341  (
343  );
344  PtrList<surfaceSymmTensorField> surfSymmTensorFlds
345  (
346  surfSymmTensorNames.size()
347  );
348  subsetSurfaceFields(subsetter, surfSymmTensorNames, surfSymmTensorFlds);
349 
350  wordList surfTensorNames(objects.names(surfaceTensorField::typeName));
351  PtrList<surfaceTensorField> surfTensorFlds(surfTensorNames.size());
352  subsetSurfaceFields(subsetter, surfTensorNames, surfTensorFlds);
353 
354 
355  // Read point fields and subset
356  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
357 
358  const pointMesh& pMesh = pointMesh::New(mesh);
359 
360  wordList pointScalarNames(objects.names(pointScalarField::typeName));
361  PtrList<pointScalarField> pointScalarFlds(pointScalarNames.size());
362  subsetPointFields(subsetter, pMesh, pointScalarNames, pointScalarFlds);
363 
364  wordList pointVectorNames(objects.names(pointVectorField::typeName));
365  PtrList<pointVectorField> pointVectorFlds(pointVectorNames.size());
366  subsetPointFields(subsetter, pMesh, pointVectorNames, pointVectorFlds);
367 
368  wordList pointSphericalTensorNames
369  (
371  );
372  PtrList<pointSphericalTensorField> pointSphericalTensorFlds
373  (
374  pointSphericalTensorNames.size()
375  );
376  subsetPointFields
377  (
378  subsetter,
379  pMesh,
380  pointSphericalTensorNames,
381  pointSphericalTensorFlds
382  );
383 
384  wordList pointSymmTensorNames
385  (
387  );
388  PtrList<pointSymmTensorField> pointSymmTensorFlds
389  (
390  pointSymmTensorNames.size()
391  );
392  subsetPointFields
393  (
394  subsetter,
395  pMesh,
396  pointSymmTensorNames,
397  pointSymmTensorFlds
398  );
399 
400  wordList pointTensorNames(objects.names(pointTensorField::typeName));
401  PtrList<pointTensorField> pointTensorFlds(pointTensorNames.size());
402  subsetPointFields(subsetter, pMesh, pointTensorNames, pointTensorFlds);
403 
404 
405  // Read dimensioned fields and subset
406  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
407 
408  typedef volScalarField::Internal dimScalType;
409  wordList scalarDimNames(objects.names(dimScalType::typeName));
410  PtrList<dimScalType> scalarDimFlds(scalarDimNames.size());
411  subsetDimensionedFields(subsetter, scalarDimNames, scalarDimFlds);
412 
413  typedef volVectorField::Internal dimVecType;
414  wordList vectorDimNames(objects.names(dimVecType::typeName));
415  PtrList<dimVecType> vectorDimFlds(vectorDimNames.size());
416  subsetDimensionedFields(subsetter, vectorDimNames, vectorDimFlds);
417 
418  typedef volSphericalTensorField::Internal dimSphereType;
419  wordList sphericalTensorDimNames
420  (
421  objects.names(dimSphereType::typeName)
422  );
423  PtrList<dimSphereType> sphericalTensorDimFlds
424  (
425  sphericalTensorDimNames.size()
426  );
427  subsetDimensionedFields
428  (
429  subsetter,
430  sphericalTensorDimNames,
431  sphericalTensorDimFlds
432  );
433 
434  typedef volSymmTensorField::Internal dimSymmTensorType;
435  wordList symmTensorDimNames(objects.names(dimSymmTensorType::typeName));
436  PtrList<dimSymmTensorType> symmTensorDimFlds(symmTensorDimNames.size());
437  subsetDimensionedFields
438  (
439  subsetter,
440  symmTensorDimNames,
441  symmTensorDimFlds
442  );
443 
444  typedef volTensorField::Internal dimTensorType;
445  wordList tensorDimNames(objects.names(dimTensorType::typeName));
446  PtrList<dimTensorType> tensorDimFlds(tensorDimNames.size());
447  subsetDimensionedFields(subsetter, tensorDimNames, tensorDimFlds);
448 
449 
450  // Write mesh and fields to new time
451  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
452 
453  if (overwrite || specifiedInstance)
454  {
455  runTime.setTime(instant(fieldsInstance), 0);
456  subsetter.subMesh().setInstance(meshInstance);
457  }
458  else
459  {
460  runTime++;
461  }
462 
463  Info<< "Writing subsetted mesh and fields to time "
464  << runTime.name() << endl;
465 
466  subsetter.subMesh().write();
467 
468  // Subsetting adds 'subset' prefix. Rename field to be like original.
469  forAll(scalarFlds, i)
470  {
471  scalarFlds[i].rename(scalarNames[i]);
472  scalarFlds[i].write();
473  }
474  forAll(vectorFlds, i)
475  {
476  vectorFlds[i].rename(vectorNames[i]);
477  vectorFlds[i].write();
478  }
479  forAll(sphericalTensorFlds, i)
480  {
481  sphericalTensorFlds[i].rename(sphericalTensorNames[i]);
482  sphericalTensorFlds[i].write();
483  }
484  forAll(symmTensorFlds, i)
485  {
486  symmTensorFlds[i].rename(symmTensorNames[i]);
487  symmTensorFlds[i].write();
488  }
489  forAll(tensorFlds, i)
490  {
491  tensorFlds[i].rename(tensorNames[i]);
492  tensorFlds[i].write();
493  }
494 
495  // Surface ones.
496  forAll(surfScalarFlds, i)
497  {
498  surfScalarFlds[i].rename(surfScalarNames[i]);
499  surfScalarFlds[i].write();
500  }
501  forAll(surfVectorFlds, i)
502  {
503  surfVectorFlds[i].rename(surfVectorNames[i]);
504  surfVectorFlds[i].write();
505  }
506  forAll(surfSphericalTensorFlds, i)
507  {
508  surfSphericalTensorFlds[i].rename(surfSphericalTensorNames[i]);
509  surfSphericalTensorFlds[i].write();
510  }
511  forAll(surfSymmTensorFlds, i)
512  {
513  surfSymmTensorFlds[i].rename(surfSymmTensorNames[i]);
514  surfSymmTensorFlds[i].write();
515  }
516  forAll(surfTensorNames, i)
517  {
518  surfTensorFlds[i].rename(surfTensorNames[i]);
519  surfTensorFlds[i].write();
520  }
521 
522  // Point ones
523  forAll(pointScalarFlds, i)
524  {
525  pointScalarFlds[i].rename(pointScalarNames[i]);
526  pointScalarFlds[i].write();
527  }
528  forAll(pointVectorFlds, i)
529  {
530  pointVectorFlds[i].rename(pointVectorNames[i]);
531  pointVectorFlds[i].write();
532  }
533  forAll(pointSphericalTensorFlds, i)
534  {
535  pointSphericalTensorFlds[i].rename(pointSphericalTensorNames[i]);
536  pointSphericalTensorFlds[i].write();
537  }
538  forAll(pointSymmTensorFlds, i)
539  {
540  pointSymmTensorFlds[i].rename(pointSymmTensorNames[i]);
541  pointSymmTensorFlds[i].write();
542  }
543  forAll(pointTensorNames, i)
544  {
545  pointTensorFlds[i].rename(pointTensorNames[i]);
546  pointTensorFlds[i].write();
547  }
548 
549  // DimensionedFields
550  forAll(scalarDimFlds, i)
551  {
552  scalarDimFlds[i].rename(scalarDimNames[i]);
553  scalarDimFlds[i].write();
554  }
555  forAll(vectorDimFlds, i)
556  {
557  vectorDimFlds[i].rename(vectorDimNames[i]);
558  vectorDimFlds[i].write();
559  }
560  forAll(sphericalTensorDimFlds, i)
561  {
562  sphericalTensorDimFlds[i].rename(sphericalTensorDimNames[i]);
563  sphericalTensorDimFlds[i].write();
564  }
565  forAll(symmTensorDimFlds, i)
566  {
567  symmTensorDimFlds[i].rename(symmTensorDimNames[i]);
568  symmTensorDimFlds[i].write();
569  }
570  forAll(tensorDimFlds, i)
571  {
572  tensorDimFlds[i].rename(tensorDimNames[i]);
573  tensorDimFlds[i].write();
574  }
575  }
576  else
577  {
578  // Write mesh to new time
579  // ~~~~~~~~~~~~~~~~~~~~~~
580 
581  if (overwrite || specifiedInstance)
582  {
583  runTime.setTime(instant(fieldsInstance), 0);
584  subsetter.subMesh().setInstance(meshInstance);
585  }
586  else
587  {
588  runTime++;
589  }
590 
591  Info<< "Writing subsetted mesh to time "
592  << runTime.name() << endl;
593 
594  subsetter.subMesh().write();
595  }
596 
597  Info<< "End\n" << endl;
598 
599  return 0;
600 }
601 
602 
603 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static const char *const typeName
Definition: Field.H:105
Generic GeometricField class.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
A collection of cell labels.
Definition: cellSet.H:51
const word & name() const
Return const reference to name.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:222
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:896
const fvMesh & subMesh() const
Return reference to subset mesh.
static tmp< VolField< Type > > interpolate(const VolField< Type > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:402
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1736
An instant of time. Contains the time value and name.
Definition: instant.H:67
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static List< word > fieldNames
Definition: globalFoam.H:46
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
error FatalError
static const char nl
Definition: Ostream.H:260
objects
Foam::argList args(argc, argv)