33 void Foam::CFCFaceToCellStencil::calcFaceBoundaryData
47 label facei = pp.start();
57 globFaces.
setSize(cFaces.size()-1);
62 if (cFaces[j] != facei)
73 else if (isA<emptyPolyPatch>(pp))
91 void Foam::CFCFaceToCellStencil::calcCellStencil
96 const label nBnd = mesh().nFaces()-mesh().nInternalFaces();
97 const labelList& own = mesh().faceOwner();
98 const labelList& nei = mesh().faceNeighbour();
105 calcFaceBoundaryData(neiGlobal);
110 boolList validBFace(mesh().nFaces()-mesh().nInternalFaces(),
true);
112 const polyBoundaryMesh&
patches = mesh().boundaryMesh();
117 if (isA<emptyPolyPatch>(pp))
119 label bFacei = pp.start()-mesh().nInternalFaces();
122 validBFace[bFacei++] =
false;
131 DynamicList<label> allGlobalFaces(100);
133 globalCellFaces.setSize(mesh().nCells());
134 forAll(globalCellFaces, celli)
136 const cell& cFaces = mesh().cells()[celli];
138 allGlobalFaces.clear();
143 label facei = cFaces[i];
147 mesh().isInternalFace(facei)
148 || validBFace[facei-mesh().nInternalFaces()]
151 allGlobalFaces.append(globalNumbering().toGlobal(facei));
158 label facei = cFaces[i];
160 if (mesh().isInternalFace(facei))
162 label nbrCelli = own[facei];
163 if (nbrCelli == celli)
165 nbrCelli = nei[facei];
167 const cell& nbrFaces = mesh().cells()[nbrCelli];
171 label nbrFacei = nbrFaces[j];
175 mesh().isInternalFace(nbrFacei)
176 || validBFace[nbrFacei-mesh().nInternalFaces()]
179 label nbrGlobalI = globalNumbering().toGlobal(nbrFacei);
182 if (
findIndex(allGlobalFaces, nbrGlobalI) == -1)
184 allGlobalFaces.append(nbrGlobalI);
192 neiGlobal[facei-mesh().nInternalFaces()];
196 label nbrGlobalI = nbrGlobalFaces[j];
199 if (
findIndex(allGlobalFaces, nbrGlobalI) == -1)
201 allGlobalFaces.append(nbrGlobalI);
207 globalCellFaces[celli] = allGlobalFaces;
230 calcCellStencil(*
this);
#define forAll(list, i)
Loop across all elements in list.
CFCFaceToCellStencil(const polyMesh &)
Construct from mesh.
void setSize(const label)
Reset size of List.
baseclass for extended cell centred addressing. Contains per cell a list of neighbouring faces in glo...
const globalIndex & globalNumbering() const
Global numbering for faces.
const polyMesh & mesh() const
label toGlobal(const label i) const
From local to global.
Mesh consisting of general polyhedral cells.
virtual const labelList & faceOwner() const
Return face owner.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
label nInternalFaces() const
const cellList & cells() const
const fvPatchList & patches
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
List< bool > boolList
Bool container classes.
List< labelList > labelListList
A List of labelList.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.