cellDecomposer.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) 2024 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "cellDecomposer.H"
30 #include "polyTopoChange.H"
31 #include "mapPolyMesh.H"
32 #include "tetDecomposer.H"
33 #include "syncTools.H"
34 #include "dummyTransform.H"
35 #include "ReadFields.H"
36 #include "surfaceFields.H"
37 #include "PackedBoolList.H"
38 #include "fvMeshTools.H"
39 #include "cellSetOption.H"
40 #include "cellBitSet.H"
41 #include "cellSet.H"
42 #include "hexMatcher.H"
43 #include "prismMatcher.H"
44 #include "pyrMatcher.H"
45 #include "tetMatcher.H"
46 
47 #include "OBJstream.H"
48 
49 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace functionObjects
54 {
55  defineTypeNameAndDebug(cellDecomposer, 0);
56  addToRunTimeSelectionTable(functionObject, cellDecomposer, dictionary);
57 }
58 }
59 
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62 
63 void Foam::functionObjects::cellDecomposer::makeMesh
64 (
65  const dictionary& dict,
66  const word& regionName
67 )
68 {
69  Info<< name() << ':' << nl
70  << " Decomposing cells to region " << regionName << endl;
71 
72  Info<< incrIndent;
73 
74  const fv::cellSetOption::selectionModeType selectionMode
75  (
77  (
78  "selectionMode",
79  dict
80  )
81  );
82 
83 
84  // Start off from existing mesh
85  polyTopoChange meshMod(mesh_);
86 
87  tetDecompPtr_.reset(new Foam::tetDecomposer(mesh_));
88 
89  // Default values
90  bitSet decomposeCell(mesh_.nCells());
91  autoPtr<bitSet> decomposeFacePtr;
93  (
95  );
96 
97 
98  switch (selectionMode)
99  {
101  {
102  Info<< indent << "- selecting all cells" << endl;
103 
104  decomposeCell = true;
105  break;
106  }
108  {
109  Info<< indent << "- selecting cells geometrically" << endl;
110 
111  decomposeCell =
112  cellBitSet::select(mesh_, dict.subDict("selection"), true);
113  break;
114  }
116  {
117  const word selectionName(dict.get<word>("cellSet"));
118 
119  Info<< indent
120  << "- selecting cells using cellSet "
121  << selectionName << endl;
122 
123  decomposeCell.set
124  (
125  cellSet(mesh_, selectionName, IOobject::MUST_READ).toc()
126  );
127  break;
128  }
130  {
131  wordRes selectionNames;
132  if
133  (
134  !dict.readIfPresent("cellZones", selectionNames)
135  || selectionNames.empty()
136  )
137  {
138  selectionNames.resize(1);
139  dict.readEntry("cellZone", selectionNames.first());
140  }
141 
142  Info<< indent
143  << "- selecting cells using cellZones "
144  << flatOutput(selectionNames) << nl;
145 
146  const auto& zones = mesh_.cellZones();
147 
148  // Also handles groups, multiple zones etc ...
149  labelList zoneIDs = zones.indices(selectionNames);
150 
151  if (zoneIDs.empty())
152  {
154  << "No matching cellZones: "
155  << flatOutput(selectionNames) << nl
156  << "Valid zones : "
157  << flatOutput(zones.names()) << nl
158  << "Valid groups: "
159  << flatOutput(zones.groupNames())
160  << nl
161  << exit(FatalError);
162  }
163 
164  if (zoneIDs.size() == 1)
165  {
166  decomposeCell.set(zones[zoneIDs.first()]);
167  // TBD: Foam::sort(cells_);
168  }
169  else
170  {
171  decomposeCell.set(zones.selection(zoneIDs).sortedToc());
172  }
173  break;
174  }
175  default:
176  {
178  << "Unsupported selectionMode "
180  << ". Valid selectionMode types are "
182  << exit(FatalError);
183  }
184  }
185 
186  word decompTypeName;
187  if (dict.readIfPresent("decomposeType", decompTypeName))
188  {
189  if (decompTypeName == "polyhedral")
190  {
191  // Automatic selection to generate hex/prism/tet only
192  // - subset decomposeCell to exclude hex/prism/tet
193  // - set faces to FACE_DIAG_QUADS
194  // Usually start with cellSetOption::smAll
195  decomposeFacePtr.reset(new bitSet(mesh_.nFaces()));
196  auto& decomposeFace = decomposeFacePtr();
197 
198  decompType = tetDecomposer::FACE_DIAG_QUADS;
199  bitSet oldDecomposeCell(decomposeCell);
200  decomposeCell = false;
201 
202  // Construct shape recognizers
203  prismMatcher prism;
204 
205  for (const label celli : oldDecomposeCell)
206  {
207  if
208  (
209  !hexMatcher::test(mesh_, celli)
210  && !tetMatcher::test(mesh_, celli)
211  && !pyrMatcher::test(mesh_, celli)
212  && !prism.isA(mesh_, celli)
213  )
214  {
215  decomposeCell.set(celli);
216  decomposeFace.set(mesh_.cells()[celli]);
217  }
218  }
219 
220  // Sync boundary info
222  (
223  mesh_,
224  decomposeFace,
225  orEqOp<unsigned int>()
226  );
227  }
228  else
229  {
230  decompType = tetDecomposer::decompositionTypeNames[decompTypeName];
231  }
232  }
233 
234 
235  // Insert all changes to create tets
236  if (decomposeFacePtr)
237  {
238  if (debug)
239  {
240  OBJstream os(mesh_.time().path()/"orig_faces.obj");
241  os.write
242  (
243  UIndirectList<face>
244  (
245  mesh_.faces(),
246  decomposeFacePtr().sortedToc()
247  )(),
248  mesh_.points(),
249  false
250  );
251  Pout<< "Written " << meshMod.faces().size()
252  << " faces to " << os.name() << endl;
253  }
254 
255  tetDecompPtr_().setRefinement
256  (
257  decompType, //Foam::tetDecomposer::FACE_CENTRE_TRIS,
258  decomposeCell,
259  decomposeFacePtr(),
260  meshMod
261  );
262  }
263  else
264  {
265  tetDecompPtr_().setRefinement
266  (
267  decompType, //Foam::tetDecomposer::FACE_CENTRE_TRIS,
268  decomposeCell,
269  meshMod
270  );
271  }
272 
273 
274  if (debug)
275  {
276  OBJstream os(mesh_.time().path()/"faces.obj");
277  os.write(meshMod.faces(), pointField(meshMod.points()), false);
278  Pout<< "Written " << meshMod.faces().size()
279  << " faces to " << os.name() << endl;
280  }
281 
282 
283  autoPtr<fvMesh> tetMeshPtr;
284 
285  mapPtr_ = meshMod.makeMesh
286  (
287  tetMeshPtr,
288  IOobject
289  (
290  regionName,
291  mesh_.facesInstance(), // same instance as original mesh
292  mesh_.time(), //? why same database? TBD.
293  IOobject::READ_IF_PRESENT, // Read fv* if present
296  ),
297  mesh_
298  );
299 
300 
301  //- Update numbering on tet-decomposition engine
302  tetDecompPtr_().updateMesh(mapPtr_());
303 
304  Info<< indent << "Writing decomposed mesh to "
305  << tetMeshPtr().objectRegistry::objectRelPath()
306  << endl;
307  tetMeshPtr().write();
308 
310 
311  // Store new mesh on object registry
312  tetMeshPtr.ptr()->polyMesh::store();
313 
314  Info<< decrIndent;
315 }
317 
318 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
319 
321 (
322  const word& name,
323  const Time& runTime,
324  const dictionary& dict
325 )
326 :
328 {
329  read(dict);
330 }
331 
332 
333 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
334 
336 {
338  {
339  // Generate new tet equivalent mesh. TBD: why meshObject?
340  //meshObjects::tetDecomposition::New(mesh, dict);
341 
342  dict_ = dict.optionalSubDict(typeName + "Coeffs");
343  dict_.readEntry("mapRegion", mapRegion_);
344  dict_.readEntry("fields", fieldNames_);
345  makeMesh(dict_, mapRegion_);
346  }
348  return true;
349 }
350 
351 
353 {
354  Log << type() << " " << name() << " execute:" << nl;
355 
356  bool ok = false;
357 
358  if (mesh_.changing())
359  {
360  // Re-do mesh. Currently rebuilds whole mesh. Could move points only
361  // for mesh motion.
362  tetDecompPtr_.clear();
363  mapPtr_.clear();
364  const_cast<Time&>(this->mesh_.time()).erase(mapRegion_);
365  makeMesh(dict_, mapRegion_);
366  }
367 
368 
369  // Look up
370  ok = mapFieldType<scalar>() || ok;
371  ok = mapFieldType<vector>() || ok;
372  ok = mapFieldType<sphericalTensor>() || ok;
373  ok = mapFieldType<symmTensor>() || ok;
374  ok = mapFieldType<tensor>() || ok;
375 
376  if (log)
377  {
378  if (!ok)
379  {
380  Info<< " none" << nl;
381  }
382 
383  Info<< endl;
384  }
385  return true;
386 }
387 
388 
390 {
391  Log << type() << " " << name() << " write:" << nl;
392 
393  bool ok = false;
394 
395  ok = writeFieldType<scalar>() || ok;
396  ok = writeFieldType<vector>() || ok;
397  ok = writeFieldType<sphericalTensor>() || ok;
398  ok = writeFieldType<symmTensor>() || ok;
399  ok = writeFieldType<tensor>() || ok;
400 
401  if (log)
402  {
403  if (!ok)
404  {
405  Info<< " none" << nl;
406  }
407 
408  Info<< endl;
409  }
410 
411  return true;
412 }
413 
414 
415 // ************************************************************************* //
Foam::surfaceFields.
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition: syncTools.H:465
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
const labelIOList & zoneIDs
Definition: correctPhi.H:59
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
dimensionedScalar log(const dimensionedScalar &ds)
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:859
srcOptions erase("case")
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const cellList & cells() const
Field reading functions for post-processing utilities.
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for HEX. (6 quad)
Definition: hexMatcher.C:80
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
label nFaces() const noexcept
Number of mesh faces.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
static const Enum< selectionModeType > selectionModeTypeNames_
List of selection mode type names.
Decomposes polyMesh into tets (or pyramids)
Definition: tetDecomposer.H:60
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
static const Enum< decompositionType > decompositionTypeNames
Definition: tetDecomposer.H:74
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
Dummy transform to be used with syncTools.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
cellDecomposer(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Reading is optional [identical to LAZY_READ].
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
Definition: fvMeshTools.C:1212
static bitSet select(const polyMesh &mesh, const dictionary &dict, const bool verbosity=false)
Return a cell selection according to the dictionary specification of actions.
Definition: cellBitSet.C:84
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
selectionModeType
Enumeration for selection mode types.
virtual bool read(const dictionary &dict)
Read the cellDecomposer data.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
label nCells() const noexcept
Number of mesh cells.
#define Log
Definition: PDRblock.C:28
Automatically write from objectRegistry::writeObject()
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
List< label > labelList
A List of labels.
Definition: List.H:62
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for PYR. (4 tri, 1 quad)
Definition: pyrMatcher.C:107
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for TET. (4 tri)
Definition: tetMatcher.C:82
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
Request registration (bool: true)
List< label > toc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:158
const fvMesh & mesh_
Reference to the fvMesh.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
static int debug
Flag to execute debug content.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225