CFCFaceToCellStencil.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) 2011-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 "CFCFaceToCellStencil.H"
27 #include "syncTools.H"
28 #include "emptyPolyPatch.H"
29 #include "dummyTransform.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::CFCFaceToCellStencil::calcFaceBoundaryData
34 (
35  labelListList& neiGlobal
36 ) const
37 {
38  const polyBoundaryMesh& patches = mesh().boundaryMesh();
39  const label nBnd = mesh().nFaces()-mesh().nInternalFaces();
40  const labelList& own = mesh().faceOwner();
41 
42  neiGlobal.setSize(nBnd);
43 
44  forAll(patches, patchi)
45  {
46  const polyPatch& pp = patches[patchi];
47  label facei = pp.start();
48 
49  if (pp.coupled())
50  {
51  // For coupled faces get the faces of the cell on the other side
52  forAll(pp, i)
53  {
54  const labelList& cFaces = mesh().cells()[own[facei]];
55 
56  labelList& globFaces = neiGlobal[facei-mesh().nInternalFaces()];
57  globFaces.setSize(cFaces.size()-1);
58  label globI = 0;
59 
60  forAll(cFaces, j)
61  {
62  if (cFaces[j] != facei)
63  {
64  globFaces[globI++] = globalNumbering().toGlobal
65  (
66  cFaces[j]
67  );
68  }
69  }
70  facei++;
71  }
72  }
73  else if (isA<emptyPolyPatch>(pp))
74  {}
75  else
76  {
77  // Do nothing since face itself already in stencil
78  }
79  }
80 
82  (
83  mesh(),
84  neiGlobal,
85  eqOp<labelList>(),
86  dummyTransform()
87  );
88 }
89 
90 
91 void Foam::CFCFaceToCellStencil::calcCellStencil
92 (
93  labelListList& globalCellFaces
94 ) const
95 {
96  const label nBnd = mesh().nFaces()-mesh().nInternalFaces();
97  const labelList& own = mesh().faceOwner();
98  const labelList& nei = mesh().faceNeighbour();
99 
100 
101  // Calculate faces of coupled neighbour (in global numbering)
102  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
103 
104  labelListList neiGlobal(nBnd);
105  calcFaceBoundaryData(neiGlobal);
106 
107 
108 
109  // Non-empty boundary faces
110  boolList validBFace(mesh().nFaces()-mesh().nInternalFaces(), true);
111 
112  const polyBoundaryMesh& patches = mesh().boundaryMesh();
113  forAll(patches, patchi)
114  {
115  const polyPatch& pp = patches[patchi];
116 
117  if (isA<emptyPolyPatch>(pp))
118  {
119  label bFacei = pp.start()-mesh().nInternalFaces();
120  forAll(pp, i)
121  {
122  validBFace[bFacei++] = false;
123  }
124  }
125  }
126 
127 
128  // Determine faces of cellCells in global numbering
129  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
130 
131  DynamicList<label> allGlobalFaces(100);
132 
133  globalCellFaces.setSize(mesh().nCells());
134  forAll(globalCellFaces, celli)
135  {
136  const cell& cFaces = mesh().cells()[celli];
137 
138  allGlobalFaces.clear();
139 
140  // My faces first
141  forAll(cFaces, i)
142  {
143  label facei = cFaces[i];
144 
145  if
146  (
147  mesh().isInternalFace(facei)
148  || validBFace[facei-mesh().nInternalFaces()]
149  )
150  {
151  allGlobalFaces.append(globalNumbering().toGlobal(facei));
152  }
153  }
154 
155  // faces of neighbouring cells second
156  forAll(cFaces, i)
157  {
158  label facei = cFaces[i];
159 
160  if (mesh().isInternalFace(facei))
161  {
162  label nbrCelli = own[facei];
163  if (nbrCelli == celli)
164  {
165  nbrCelli = nei[facei];
166  }
167  const cell& nbrFaces = mesh().cells()[nbrCelli];
168 
169  forAll(nbrFaces, j)
170  {
171  label nbrFacei = nbrFaces[j];
172 
173  if
174  (
175  mesh().isInternalFace(nbrFacei)
176  || validBFace[nbrFacei-mesh().nInternalFaces()]
177  )
178  {
179  label nbrGlobalI = globalNumbering().toGlobal(nbrFacei);
180 
181  // Check if already there. Note:should use hashset?
182  if (findIndex(allGlobalFaces, nbrGlobalI) == -1)
183  {
184  allGlobalFaces.append(nbrGlobalI);
185  }
186  }
187  }
188  }
189  else
190  {
191  const labelList& nbrGlobalFaces =
192  neiGlobal[facei-mesh().nInternalFaces()];
193 
194  forAll(nbrGlobalFaces, j)
195  {
196  label nbrGlobalI = nbrGlobalFaces[j];
197 
198  // Check if already there. Note:should use hashset?
199  if (findIndex(allGlobalFaces, nbrGlobalI) == -1)
200  {
201  allGlobalFaces.append(nbrGlobalI);
202  }
203  }
204  }
205  }
206 
207  globalCellFaces[celli] = allGlobalFaces;
208  // Pout<< "** cell:" << celli
209  // << " at:" << mesh().cellCentres()[celli]
210  // << endl;
211  // const labelList& globalFaces = globalCellFaces[celli];
212  // forAll(globalFaces, i)
213  //{
214  // label facei = globalNumbering().toLocal(globalFaces[i]);
215  // Pout<< " face:" << facei
216  // << " at:" << mesh().faceCentres()[facei]
217  // << endl;
218  //}
219  }
220 }
221 
222 
223 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
224 
226 :
227  faceToCellStencil(mesh)
228 {
229  // Calculate per cell the (face) connected cells (in global numbering)
230  calcCellStencil(*this);
231 }
232 
233 
234 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
static void syncBoundaryFaceList(const polyMesh &, UList< T > &, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
CFCFaceToCellStencil(const polyMesh &)
Construct from mesh.
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
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1047
label nFaces() const
const globalIndex & globalNumbering() const
Global numbering for faces.
baseclass for extended cell centred addressing. Contains per cell a list of neighbouring faces in glo...
const cellList & cells() const
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
Dummy transform to be used with syncTools.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:124
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1041
List< label > labelList
A List of labels.
Definition: labelList.H:56
const polyMesh & mesh() const
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void setSize(const label)
Reset size of List.
Definition: List.C:281
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
label patchi
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74