structuredRenumber.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) 2012-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "structuredRenumber.H"
27 #include "fvMeshSubset.H"
28 #include "OppositeFaceCellWave.H"
29 #include "globalIndex.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 
39  (
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 Foam::structuredRenumber::structuredRenumber
50 (
51  const dictionary& renumberDict,
52  const dictionary& methodDict
53 )
54 :
55  renumberMethod(renumberDict),
56  patches_(methodDict.lookup("patches")),
57  nLayers_(methodDict.lookupOrDefault<label>("nLayers", labelMax)),
58  depthFirst_(methodDict.lookup("depthFirst")),
59  method_(renumberMethod::New(methodDict)),
60  reverse_(methodDict.lookup("reverse"))
61 {}
62 
63 
64 Foam::structuredRenumber::structuredRenumber
65 (
66  const dictionary& renumberDict
67 )
68 :
70  (
71  renumberDict,
72  renumberDict.optionalSubDict(typeName + "Coeffs")
73  )
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
79 bool Foam::structuredRenumber::layerLess::operator()
80 (
81  const label a,
82  const label b
83 )
84 {
85  const topoDistanceData& ta = distance_[a];
86  const topoDistanceData& tb = distance_[b];
87 
88  int dummy;
89 
90  if (ta.valid(dummy))
91  {
92  if (tb.valid(dummy))
93  {
94  if (depthFirst_)
95  {
96  if (ta.data() < tb.data())
97  {
98  // Sort column first
99  return true;
100  }
101  else if (ta.data() > tb.data())
102  {
103  return false;
104  }
105  else
106  {
107  // Same column. Sort according to layer
108  return ta.distance() < tb.distance();
109  }
110  }
111  else
112  {
113  if (ta.distance() < tb.distance())
114  {
115  return true;
116  }
117  else if (ta.distance() > tb.distance())
118  {
119  return false;
120  }
121  else
122  {
123  // Same layer; sort according to current values
124  return ta.data() < tb.data();
125  }
126  }
127  }
128  else
129  {
130  return true;
131  }
132  }
133  else
134  {
135  if (tb.valid(dummy))
136  {
137  return false;
138  }
139  else
140  {
141  // Both not valid; fall back to cell indices for sorting
142  return order_[a] < order_[b];
143  }
144  }
145 }
146 
147 
149 (
150  const polyMesh& mesh,
151  const pointField& points
152 ) const
153 {
154  if (points.size() != mesh.nCells())
155  {
157  << "Number of points " << points.size()
158  << " should equal the number of cells " << mesh.nCells()
159  << exit(FatalError);
160  }
161 
162  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
163  const labelHashSet patchIDs(pbm.patchSet(patches_));
164 
165  label nFaces = 0;
166  forAllConstIter(labelHashSet, patchIDs, iter)
167  {
168  nFaces += pbm[iter.key()].size();
169  }
170 
171 
172  // Extract a submesh.
173  labelHashSet patchCells(2*nFaces);
174  forAllConstIter(labelHashSet, patchIDs, iter)
175  {
176  const labelUList& fc = pbm[iter.key()].faceCells();
177  forAll(fc, i)
178  {
179  patchCells.insert(fc[i]);
180  }
181  }
182 
183  label nTotalSeeds = returnReduce(patchCells.size(), sumOp<label>());
184 
185  label nTotalCells = mesh.globalData().nTotalCells();
186  const label nLayers = nTotalCells/nTotalSeeds;
187 
188  Info<< type() << " : seeding " << nTotalSeeds
189  << " cells on (estimated) " << nLayers << " layers" << nl
190  << endl;
191 
192 
193  // Work array. Used here to temporarily store the original-to-ordered
194  // index. Later on used to store the ordered-to-original.
195  labelList orderedToOld(mesh.nCells(), -1);
196 
197  // Subset the layer of cells next to the patch
198  {
199  fvMeshSubset subsetter(dynamic_cast<const fvMesh&>(mesh));
200  subsetter.setLargeCellSubset(patchCells);
201  const fvMesh& subMesh = subsetter.subMesh();
202 
203  pointField subPoints(points, subsetter.cellMap());
204 
205  // Locally renumber the layer of cells
206  labelList subOrder(method_().renumber(subMesh, subPoints));
207 
208  labelList subOrigToOrdered(invert(subOrder.size(), subOrder));
209 
210  globalIndex globalSubCells(subOrder.size());
211 
212  // Transfer to final decomposition and convert into global numbering
213  forAll(subOrder, i)
214  {
215  orderedToOld[subsetter.cellMap()[i]] =
216  globalSubCells.toGlobal(subOrigToOrdered[i]);
217  }
218  }
219 
220 
221  // Walk sub-ordering (=column index) out.
222  labelList patchFaces(nFaces);
223  List<topoDistanceData> patchData(nFaces);
224  nFaces = 0;
225  forAllConstIter(labelHashSet, patchIDs, iter)
226  {
227  const polyPatch& pp = pbm[iter.key()];
228  const labelUList& fc = pp.faceCells();
229  forAll(fc, i)
230  {
231  patchFaces[nFaces] = pp.start()+i;
232  patchData[nFaces] = topoDistanceData
233  (
234  orderedToOld[fc[i]],// passive data: global column
235  0 // distance: layer
236  );
237  nFaces++;
238  }
239  }
240 
241  // Field on cells and faces.
242  List<topoDistanceData> cellData(mesh.nCells());
244 
245  // Propagate information inwards
247  (
248  mesh,
249  patchFaces,
250  patchData,
251  faceData,
252  cellData,
253  0
254  );
255 
256  deltaCalc.iterate(nLayers_);
257 
258  Info<< type() << " : did not visit "
259  << deltaCalc.getUnsetCells()
260  << " cells out of " << nTotalCells
261  << "; using " << method_().type() << " renumbering for these" << endl;
262 
263  // Get cell order using the method(). These values will get overwritten
264  // by any visited cell so are used only if the number of nLayers is limited.
265  labelList oldToOrdered
266  (
267  invert
268  (
269  mesh.nCells(),
270  method_().renumber(mesh, points)
271  )
272  );
273 
274  // Use specialised sorting to sorted either layers or columns first
275  // Done so that at no point we need to combine both into a single
276  // index and we might run out of label size.
278  (
279  cellData,
280  orderedToOld,
281  layerLess(depthFirst_, oldToOrdered, cellData)
282  );
283 
284  // Return furthest away cell first
285  if (reverse_)
286  {
287  reverse(orderedToOld);
288  }
289 
290  return orderedToOld;
291 }
292 
293 
294 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
label getUnsetCells() const
Get number of unvisited cells, i.e. cells that were not (yet)
Definition: FaceCellWave.C:778
virtual label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
Definition: FaceCellWave.C:951
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Version of FaceCellWave that walks through prismatic cells only.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base class for classes which manage incomplete sets of face data.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
const labelList & cellMap() const
Return cell map.
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:917
const fvMesh & subMesh() const
Return reference to subset mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
label nTotalCells() const
Return total number of cells in decomposed mesh.
Foam::polyBoundaryMesh.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the patch set corresponding to the given names.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:313
label nCells() const
label nFaces() const
Abstract base class for renumbering.
Less function class that can be used for sorting according to.
Renumbering according to mesh layers. depthFirst = true: first column gets ids 0.....
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
For use with FaceCellWave. Determines topological distance to starting faces.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
volScalarField & b
Definition: createFields.H:27
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
messageStream Info
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
static const label labelMax
Definition: label.H:62
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
static const char nl
Definition: Ostream.H:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488