wallLayerCells.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) 2011-2016 OpenFOAM Foundation
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 "wallLayerCells.H"
29 #include "DynamicList.H"
30 #include "MeshWave.H"
31 #include "wallNormalInfo.H"
32 #include "OFstream.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(wallLayerCells, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 bool Foam::wallLayerCells::usesCoupledPatch(const label celli) const
45 {
46  const polyBoundaryMesh& patches = mesh().boundaryMesh();
47 
48  const cell& cFaces = mesh().cells()[celli];
49 
50  forAll(cFaces, cFacei)
51  {
52  label facei = cFaces[cFacei];
53 
54  label patchID = patches.whichPatch(facei);
55 
56  if ((patchID >= 0) && (patches[patchID].coupled()))
57  {
58  return true;
59  }
60  }
61  return false;
62 }
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
67 (
68  const polyMesh& mesh,
69  const List<word>& patchNames,
70  const label nLayers
71 )
72 :
73  edgeVertex(mesh),
74  List<refineCell>()
75 {
76  // Find out cells connected to walls.
77 
79 
80  // Make map from name to local patch ID
81  HashTable<label> patchNameToIndex(patches.size());
82 
83  forAll(patches, patchi)
84  {
85  patchNameToIndex.insert(patches[patchi].name(), patchi);
86  }
87 
88 
89  // Count size of walls to set
90  label nWalls = 0;
91 
92  forAll(patchNames, patchNameI)
93  {
94  const word& name = patchNames[patchNameI];
95 
96  if (patchNameToIndex.found(name))
97  {
98  label patchi = patchNameToIndex[name];
99 
100  nWalls += patches[patchi].size();
101  }
102  }
103 
104  // Allocate storage for start of wave on faces
105  List<wallNormalInfo> changedFacesInfo(nWalls);
106  labelList changedFaces(nWalls);
107 
108  // Fill changedFaces info
109  label nChangedFaces = 0;
110 
111  forAll(patchNames, patchNameI)
112  {
113  const word& name = patchNames[patchNameI];
114 
115  if (patchNameToIndex.found(name))
116  {
117  label patchi = patchNameToIndex[name];
118 
119  const polyPatch& pp = patches[patchi];
120 
121  forAll(pp, patchFacei)
122  {
123  label meshFacei = pp.start() + patchFacei;
124 
125  changedFaces[nChangedFaces] = meshFacei;
126 
127  // Set transported information to the wall normal.
128  const vector& norm = pp.faceNormals()[patchFacei];
129 
130  changedFacesInfo[nChangedFaces] = wallNormalInfo(norm);
131 
132  nChangedFaces++;
133  }
134  }
135  }
136 
137 
138  // Do a wave of nLayers, transporting the index in patchNames
139  // (cannot use local patchIDs since we might get info from neighbouring
140  // processor)
141 
142  MeshWave<wallNormalInfo> regionCalc
143  (
144  mesh,
145  changedFaces,
146  changedFacesInfo,
147  0
148  );
149 
150  regionCalc.iterate(nLayers);
151 
152 
153  // Now regionCalc should hold info on faces that are reachable from
154  // changedFaces within nLayers iterations. We use face info since that is
155  // guaranteed to be consistent across processor boundaries.
156 
157  const List<wallNormalInfo>& faceInfo = regionCalc.allFaceInfo();
158 
159  if (debug)
160  {
161  InfoInFunction << "Dumping selected faces to selectedFaces.obj" << endl;
162 
163  OFstream fcStream("selectedFaces.obj");
164 
165  label vertI = 0;
166 
167  forAll(faceInfo, facei)
168  {
169  const wallNormalInfo& info = faceInfo[facei];
170 
171  if (info.valid(regionCalc.data()))
172  {
173  const face& f = mesh.faces()[facei];
174 
175  point mid(Zero);
176 
177  forAll(f, fp)
178  {
179  mid += mesh.points()[f[fp]];
180  }
181  mid /= f.size();
182 
183  fcStream
184  << "v " << mid.x() << ' ' << mid.y() << ' ' << mid.z()
185  << endl;
186  vertI++;
187 
188  point end(mid + info.normal());
189 
190  fcStream
191  << "v " << end.x() << ' ' << end.y() << ' ' << end.z()
192  << endl;
193  vertI++;
194 
195  fcStream << "l " << vertI << ' ' <<vertI-1 << endl;
196  }
197  }
198  }
199 
200 
201  //
202  // Copy meshWave information to List<refineCell>
203  //
204 
205  // Estimate size
206 
207  DynamicList<refineCell> refineCells(3*nWalls);
208 
209  const List<wallNormalInfo>& cellInfo = regionCalc.allCellInfo();
210 
211  forAll(cellInfo, celli)
212  {
213  const wallNormalInfo& info = cellInfo[celli];
214 
215  if (info.valid(regionCalc.data()) && !usesCoupledPatch(celli))
216  {
217  refineCells.append(refineCell(celli, info.normal()));
218  }
219  }
220 
221  // Transfer refineCells storage to this.
222  transfer(refineCells);
223 }
224 
225 
226 // ************************************************************************* //
wallLayerCells(const polyMesh &mesh, const List< word > &patchNames, const label nLayers)
Construct from components.
void transfer(List< refineCell > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const cellList & cells() const
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
wordList patchNames(nPatches)
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
const polyMesh & mesh() const
Definition: edgeVertex.H:115
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
vector point
Point is a vector.
Definition: point.H:37
const polyBoundaryMesh & patches
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:435
List< label > labelList
A List of labels.
Definition: List.H:62
bool coupled
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127