structuredRenumber.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) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2018-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 \*---------------------------------------------------------------------------*/
28 
29 #include "structuredRenumber.H"
31 #include "topoDistanceData.H"
32 #include "fvMeshSubset.H"
33 #include "OppositeFaceCellWave.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(structuredRenumber, 0);
40 
42  (
43  renumberMethod,
44  structuredRenumber,
45  dictionary
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 Foam::structuredRenumber::structuredRenumber
53 (
54  const dictionary& dict
55 )
56 :
58  coeffsDict_(dict.optionalSubDict(typeName + "Coeffs")),
59  patches_(coeffsDict_.get<wordRes>("patches")),
60  nLayers_(coeffsDict_.getOrDefault<label>("nLayers", labelMax)),
61  depthFirst_(coeffsDict_.get<bool>("depthFirst")),
62  reverse_(coeffsDict_.getOrDefault("reverse", false)),
63  method_(renumberMethod::New(coeffsDict_))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 bool Foam::structuredRenumber::layerLess::operator()
70 (
71  const label a,
72  const label b
73 ) const
74 {
75  const auto& ta = distance_[a];
76  const auto& tb = distance_[b];
77 
78  int dummy;
79 
80  if (ta.valid(dummy))
81  {
82  if (tb.valid(dummy))
83  {
84  if (depthFirst_)
85  {
86  if (ta.data() < tb.data())
87  {
88  // Sort column first
89  return true;
90  }
91  else if (ta.data() > tb.data())
92  {
93  return false;
94  }
95  else
96  {
97  // Same column. Sort according to layer
98  return ta.distance() < tb.distance();
99  }
100  }
101  else
102  {
103  if (ta.distance() < tb.distance())
104  {
105  return true;
106  }
107  else if (ta.distance() > tb.distance())
108  {
109  return false;
110  }
111  else
112  {
113  // Same layer; sort according to current values
114  return ta.data() < tb.data();
115  }
116  }
117  }
118  else
119  {
120  return true;
121  }
122  }
123  else
124  {
125  if (tb.valid(dummy))
126  {
127  return false;
128  }
129  else
130  {
131  // Both not valid; fall back to cell indices for sorting
132  return order_[a] < order_[b];
133  }
134  }
135 }
136 
137 
139 (
140  const polyMesh& mesh
141 ) const
142 {
143  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
144  const labelHashSet patchIDs(pbm.patchSet(patches_));
145 
146  label nFaces = 0;
147  for (const label patchi : patchIDs)
148  {
149  nFaces += pbm[patchi].size();
150  }
151 
152 
153  // Extract a submesh.
154  labelHashSet patchCells(2*nFaces);
155  for (const label patchId : patchIDs)
156  {
157  patchCells.insert(pbm[patchId].faceCells());
158  }
159 
160  label nTotalSeeds = returnReduce(patchCells.size(), sumOp<label>());
161 
162  label nTotalCells = mesh.globalData().nTotalCells();
163  const label nLayers = nTotalCells/nTotalSeeds;
164 
165  Info<< type() << " : seeding " << nTotalSeeds
166  << " cells on (estimated) " << nLayers << " layers" << nl
167  << endl;
168 
169 
170  // Work array. Used here to temporarily store the original-to-ordered
171  // index. Later on used to store the ordered-to-original.
172  labelList orderedToOld(mesh.nCells(), -1);
173 
174  // Subset the layer of cells next to the patch
175  {
176  fvMeshSubset subsetter
177  (
178  dynamic_cast<const fvMesh&>(mesh),
179  patchCells
180  );
181 
182  const labelList& cellMap = subsetter.cellMap();
183 
184  // Locally renumber the layer of cells
185  labelList subOrder = method_().renumber(subsetter.subMesh());
186 
187  labelList subOrigToOrdered(invert(subOrder.size(), subOrder));
188 
189  const globalIndex globalSubCells(subOrder.size());
190 
191  // Transfer to final decomposition and convert into global numbering
192  forAll(subOrder, i)
193  {
194  orderedToOld[cellMap[i]] =
195  globalSubCells.toGlobal(subOrigToOrdered[i]);
196  }
197  }
198 
199 
200  // Walk sub-ordering (=column index) out.
201  labelList patchFaces(nFaces);
202  List<topoDistanceData<label>> patchData(nFaces);
203  nFaces = 0;
204  for (const label patchi : patchIDs)
205  {
206  const polyPatch& pp = pbm[patchi];
207  const labelUList& fc = pp.faceCells();
208  forAll(fc, i)
209  {
210  patchFaces[nFaces] = pp.start()+i;
211  patchData[nFaces] = topoDistanceData<label>
212  (
213  0, // distance: layer
214  orderedToOld[fc[i]] // passive data: global column
215  );
216  nFaces++;
217  }
218  }
219 
220  // Field on cells and faces.
221  List<topoDistanceData<label>> cellData(mesh.nCells());
222  List<topoDistanceData<label>> faceData(mesh.nFaces());
223 
224  // Propagate information inwards
225  OppositeFaceCellWave<topoDistanceData<label>> deltaCalc
226  (
227  mesh,
228  patchFaces,
229  patchData,
230  faceData,
231  cellData,
232  0
233  );
234 
235  deltaCalc.iterate(nLayers_);
236 
237  Info<< type() << " : did not visit "
238  << deltaCalc.nUnvisitedCells()
239  << " cells out of " << nTotalCells
240  << "; using " << method_().type() << " renumbering for these" << endl;
241 
242  // Get cell order using the method(). These values will get overwritten
243  // by any visited cell so are used only if the number of nLayers is limited.
244  labelList oldToOrdered
245  (
246  invert(mesh.nCells(), method_().renumber(mesh))
247  );
248 
249  // Use specialised sorting to sorted either layers or columns first
250  // Done so that at no point we need to combine both into a single
251  // index and we might run out of label size.
253  (
254  cellData,
255  orderedToOld,
256  layerLess(depthFirst_, oldToOrdered, cellData)
257  );
258 
259  // Return furthest away cell first
260  if (reverse_)
261  {
262  Foam::reverse(orderedToOld);
263  }
264 
265  return orderedToOld;
266 }
267 
268 
269 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
label patchId(-1)
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const polyBoundaryMesh & pbm
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual labelList renumber(const polyMesh &mesh) const
Return the cell visit order (from ordered back to original cell id)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
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
Abstract base class for renumbering.
dynamicFvMesh & mesh
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition: UListI.H:521
defineTypeNameAndDebug(combustionModel, 0)
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:30
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
List< label > labelList
A List of labels.
Definition: List.H:62
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)