zoneCellStencils.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) 2020 DLR
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 "zoneCellStencils.H"
29 #include "syncTools.H"
30 #include "dummyTransform.H"
31 #include "emptyPolyPatch.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(zoneCellStencils, 0);
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
45 {
47 
48  label nNonEmpty = 0;
49 
50  for (const polyPatch& pp : patches)
51  {
52  if (!isA<emptyPolyPatch>(pp))
53  {
54  nNonEmpty += pp.size();
55  }
56  }
57  labelList nonEmptyFaces(nNonEmpty);
58  nNonEmpty = 0;
59 
60  for (const polyPatch& pp : patches)
61  {
62  if (!isA<emptyPolyPatch>(pp))
63  {
64  label facei = pp.start();
65 
66  forAll(pp, i)
67  {
68  nonEmptyFaces[nNonEmpty++] = facei++;
69  }
70  }
71  }
72 
74  (
76  (
77  mesh_.faces(),
78  nonEmptyFaces
79  ),
80  mesh_.points()
81  );
82 }
83 
84 
87 {
88  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
89 
90  label nCoupled = 0;
91 
92  for (const polyPatch& pp : patches)
93  {
94  if (pp.coupled())
95  {
96  nCoupled += pp.size();
97  }
98  }
99  labelList coupledFaces(nCoupled);
100  nCoupled = 0;
101 
102  for (const polyPatch& pp : patches)
103  {
104  if (pp.coupled())
105  {
106  label facei = pp.start();
107 
108  forAll(pp, i)
109  {
110  coupledFaces[nCoupled++] = facei++;
111  }
112  }
113  }
114 
116  (
118  (
119  mesh_.faces(),
120  coupledFaces
121  ),
122  mesh_.points()
123  );
124 }
125 
126 
128 {
129  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
130 
131  isValidBFace.setSize(mesh_.nBoundaryFaces());
132 
133  isValidBFace = true;
134 
135  for (const polyPatch& pp : patches)
136  {
137  if (pp.coupled() || isA<emptyPolyPatch>(pp))
138  {
139  label bFacei = pp.start()-mesh_.nInternalFaces();
140  forAll(pp, i)
141  {
142  isValidBFace[bFacei++] = false;
143  }
144  }
145  }
146 }
147 
148 
150 (
151  const label globalI,
152  const labelList& pGlobals,
153  labelList& cCells
154 )
155 {
157  for (const label celli : cCells)
158  {
159  if (celli != globalI)
160  {
161  set.insert(celli);
162  }
163  }
164 
165  for (const label celli : pGlobals)
166  {
167  if (celli != globalI)
168  {
169  set.insert(celli);
170  }
171  }
172 
173  cCells.setSize(set.size()+1);
174  label n = 0;
175  cCells[n++] = globalI;
176 
177  for (const label seti : set)
178  {
179  cCells[n++] = seti;
180  }
181 }
182 
183 
185 (
186  const label exclude0,
187  const label exclude1,
188  const boolList& isValidBFace,
189  const labelList& faceLabels,
190  labelHashSet& globals
191 ) const
192 {
193  const labelList& own = mesh_.faceOwner();
194  const labelList& nei = mesh_.faceNeighbour();
195 
196  forAll(faceLabels, i)
197  {
198  label facei = faceLabels[i];
199 
200  label globalOwn = globalNumbering().toGlobal(own[facei]);
201  if (globalOwn != exclude0 && globalOwn != exclude1)
202  {
203  globals.insert(globalOwn);
204  }
205 
206  if (mesh_.isInternalFace(facei))
207  {
208  label globalNei = globalNumbering().toGlobal(nei[facei]);
209  if (globalNei != exclude0 && globalNei != exclude1)
210  {
211  globals.insert(globalNei);
212  }
213  }
214  else
215  {
216  const label bFacei = facei-mesh_.nInternalFaces();
217 
218  if (isValidBFace[bFacei])
219  {
220  label globalI = globalNumbering().toGlobal
221  (
222  mesh_.nCells() + bFacei
223  );
224 
225  if (globalI != exclude0 && globalI != exclude1)
226  {
227  globals.insert(globalI);
228  }
229  }
230  }
231  }
232 }
233 
234 
236 (
237  const boolList& isValidBFace,
238  const labelList& faceLabels,
239  labelHashSet& globals
240 ) const
241 {
242  globals.clear();
243 
244  insertFaceCells
245  (
246  -1,
247  -1,
248  isValidBFace,
249  faceLabels,
250  globals
251  );
252 
253  return globals.toc();
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
258 
260 :
262  labelListList(mesh.nCells()),
263  needComm_(),
264  globalNumbering_(mesh_.nCells() + mesh_.nBoundaryFaces())
265 {}
266 
267 
268 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
269 
271 {
272  if (mesh_.topoChanging())
273  {
274  globalNumbering_ =
275  globalIndex(mesh_.nCells() + mesh_.nBoundaryFaces());
276 
277  static_cast<labelListList&>(*this) = labelListList(mesh_.nCells());
278  needComm_.clear();
279  }
280 }
281 
282 
283 // ************************************************************************* //
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::zoneCellStencils::nonEmptyFacesPatch
autoPtr< indirectPrimitivePatch > nonEmptyFacesPatch() const
Return patch of all coupled faces.
Definition: zoneCellStencils.C:44
Foam::zoneCellStencils::allCoupledFacesPatch
autoPtr< indirectPrimitivePatch > allCoupledFacesPatch() const
Return patch of all coupled faces.
Definition: zoneCellStencils.C:86
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::HashTable::toc
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:121
Foam::zoneCellStencils::zoneCellStencils
zoneCellStencils(const fvMesh &)
Construct from all cells and boundary faces.
Definition: zoneCellStencils.C:259
Foam::zoneCellStencils::validBoundaryFaces
void validBoundaryFaces(boolList &isValidBFace) const
Valid boundary faces (not empty and not coupled)
Definition: zoneCellStencils.C:127
Foam::zoneCellStencils::merge
static void merge(const label globalI, const labelList &pGlobals, labelList &cCells)
Merge two lists and guarantee globalI is first.
Definition: zoneCellStencils.C:150
dummyTransform.H
Dummy transform to be used with syncTools.
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::HashSet< label, Hash< label > >
Foam::UpdateableMeshObject
Definition: MeshObject.H:241
syncTools.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::zoneCellStencils::insertFaceCells
void insertFaceCells(const label exclude0, const label exclude1, const boolList &nonEmptyFace, const labelList &faceLabels, labelHashSet &globals) const
Collect cell neighbours of faces in global numbering.
Definition: zoneCellStencils.C:185
n
label n
Definition: TABSMDCalcMethod2.H:31
zoneCellStencils.H
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:56
Foam::MeshObject< fvMesh, UpdateableMeshObject, zoneCellStencils >::mesh_
const fvMesh & mesh_
Definition: MeshObject.H:96
Foam::zoneCellStencils
base class for cell stencil in a narrow band
Definition: zoneCellStencils.H:58
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
emptyPolyPatch.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::zoneCellStencils::calcFaceCells
labelList calcFaceCells(const boolList &nonEmptyFace, const labelList &faceLabels, labelHashSet &globals) const
Collect cell neighbours of faces in global numbering.
Definition: zoneCellStencils.C:236
Foam::HashTable::clear
void clear()
Clear all entries from table.
Definition: HashTable.C:630
Foam::List< label >
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:181
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:88
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::zoneCellStencils::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Definition: zoneCellStencils.C:270