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 | Copyright (C) 2012-2015 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 "FaceCellWave.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.subDict(typeName + "Coeffs")),
56  patches_(methodDict_.lookup("patches")),
57  //nLayers_(readLabel(methodDict_.lookup("nLayers"))),
58  depthFirst_(methodDict_.lookup("depthFirst")),
59  method_(renumberMethod::New(methodDict_)),
60  reverse_(methodDict_.lookup("reverse"))
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
67 (
68  const polyMesh& mesh,
69  const pointField& points
70 ) const
71 {
72  if (points.size() != mesh.nCells())
73  {
75  << "Number of points " << points.size()
76  << " should equal the number of cells " << mesh.nCells()
77  << exit(FatalError);
78  }
79 
80  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
81  const labelHashSet patchIDs(pbm.patchSet(patches_));
82 
83  label nFaces = 0;
84  forAllConstIter(labelHashSet, patchIDs, iter)
85  {
86  nFaces += pbm[iter.key()].size();
87  }
88 
89 
90  // Extract a submesh.
91  labelHashSet patchCells(2*nFaces);
92  forAllConstIter(labelHashSet, patchIDs, iter)
93  {
94  const labelUList& fc = pbm[iter.key()].faceCells();
95  forAll(fc, i)
96  {
97  patchCells.insert(fc[i]);
98  }
99  }
100 
101  label nTotalSeeds = returnReduce(patchCells.size(), sumOp<label>());
102 
103  label nTotalCells = mesh.globalData().nTotalCells();
104  const label nLayers = nTotalCells/nTotalSeeds;
105 
106  Info<< type() << " : seeding " << nTotalSeeds
107  << " cells on " << nLayers << " layers" << nl
108  << endl;
109 
110 
111  // Avoid subsetMesh, FaceCellWave going through proc boundaries
112  bool oldParRun = Pstream::parRun();
113  Pstream::parRun() = false;
114 
115 
116  // Work array. Used here to temporarily store the original-to-ordered
117  // index. Later on used to store the ordered-to-original.
118  labelList orderedToOld(points.size(), -1);
119 
120  // Subset the layer of cells next to the patch
121  {
122  fvMeshSubset subsetter(dynamic_cast<const fvMesh&>(mesh));
123  subsetter.setLargeCellSubset(patchCells);
124  const fvMesh& subMesh = subsetter.subMesh();
125 
126  pointField subPoints(points, subsetter.cellMap());
127 
128  // Decompose the layer of cells
129  labelList subOrder(method_().renumber(subMesh, subPoints));
130 
131  labelList subOrigToOrdered(invert(subOrder.size(), subOrder));
132 
133  // Transfer to final decomposition
134  forAll(subOrder, i)
135  {
136  orderedToOld[subsetter.cellMap()[i]] = subOrigToOrdered[i];
137  }
138  }
139 
140 
141  // Walk out.
142  labelList patchFaces(nFaces);
143  List<topoDistanceData> patchData(nFaces);
144  nFaces = 0;
145  forAllConstIter(labelHashSet, patchIDs, iter)
146  {
147  const polyPatch& pp = pbm[iter.key()];
148  const labelUList& fc = pp.faceCells();
149  forAll(fc, i)
150  {
151  patchFaces[nFaces] = pp.start()+i;
152  patchData[nFaces] = topoDistanceData
153  (
154  orderedToOld[fc[i]],// passive data: order of originating face
155  0 // distance: layer
156  );
157  nFaces++;
158  }
159  }
160 
161  // Field on cells and faces.
162  List<topoDistanceData> cellData(mesh.nCells());
163  List<topoDistanceData> faceData(mesh.nFaces());
164 
165  // Propagate information inwards
167  (
168  mesh,
169  patchFaces,
170  patchData,
171  faceData,
172  cellData,
173  nTotalCells+1
174  );
175 
176 
177  Pstream::parRun() = oldParRun;
178 
179 
180  // And extract.
181  // Note that distance is distance from face so starts at 1.
182  bool haveWarned = false;
183  forAll(orderedToOld, cellI)
184  {
185  if (!cellData[cellI].valid(deltaCalc.data()))
186  {
187  if (!haveWarned)
188  {
190  << "Did not visit some cells, e.g. cell " << cellI
191  << " at " << mesh.cellCentres()[cellI] << endl
192  << "Assigning these cells to domain 0." << endl;
193  haveWarned = true;
194  }
195  orderedToOld[cellI] = 0;
196  }
197  else
198  {
199  label layerI = cellData[cellI].distance();
200  if (depthFirst_)
201  {
202  orderedToOld[nLayers*cellData[cellI].data()+layerI] = cellI;
203  }
204  else
205  {
206  orderedToOld[cellData[cellI].data()+nLayers*layerI] = cellI;
207  }
208  }
209  }
210 
211  // Return furthest away cell first
212  if (reverse_)
213  {
214  reverse(orderedToOld);
215  }
216 
217  return orderedToOld;
218 }
219 
220 
221 // ************************************************************************* //
Foam::fvMeshSubset::setLargeCellSubset
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:795
topoDistanceData.H
Foam::topoDistanceData
For use with FaceCellWave. Determines topological distance to starting faces.
Definition: topoDistanceData.H:53
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
fvMeshSubset.H
structuredRenumber.H
Foam::FaceCellWave::data
const TrackingData & data() const
Additional data to be passed into container.
Definition: FaceCellWave.H:331
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::structuredRenumber::structuredRenumber
structuredRenumber(const structuredRenumber &)
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
Foam::HashSet< label, Hash< label > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::fvMeshSubset::cellMap
const labelList & cellMap() const
Return cell map.
Definition: fvMeshSubset.C:1542
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::polyPatch::faceCells
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::renumberMethod
Abstract base class for renumbering.
Definition: renumberMethod.H:48
Foam::FaceCellWave
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:75
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
FaceCellWave.H
Foam::UList
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::fvMeshSubset::subMesh
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubset.C:1472
Foam::invert
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::structuredRenumber::renumber
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
Definition: structuredRenumber.H:98
Foam::polyBoundaryMesh::patchSet
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.
Definition: polyBoundaryMesh.C:750
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
renumber
static void renumber(const labelList &map, labelList &elems)
Definition: reconstructParMesh.C:61
Foam::reverse
void reverse(UList< T > &, const label n)
Definition: UListI.H:322