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