PDRmeshArrays.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) 2016 Shell Research Ltd.
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "PDRmeshArrays.H"
30 #include "PDRblock.H"
31 #include "polyMesh.H"
32 #include "Time.H"
33 #include "IjkField.H"
34 
35 // Notes
36 //
37 // Determines the face and cell numbers of all faces and cells in the
38 // central rectangular region where CAD_PDR operates. First,
39 // "points" is read and the coordinates (by which I mean here the
40 // indices in the x, y and z coordinate arrays) are determined. Then
41 // "faces" is read and for each the coordinates of the lower- x,y,z
42 // corner are determioned, also the orientation (X, Y or Z).
43 // (Orientation in the sense of e.g. + or -x is not noted.) The files
44 // "owner" and "neighbour" specify the six faces around each cell, so
45 // from these the coordinates of the cells are determined.
46 //
47 // Full checks are made that the mesh in the central region is consistent
48 // with CAD_PDR's mesh specified by the PDRmeshSpec file.
49 //
50 // Eventually, when writing out results, we shall work through the
51 // full list of cells, writing default values for any cells that are
52 // not in the central regtion.
53 
54 
55 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
56 
57 Foam::scalar Foam::PDRmeshArrays::gridPointRelTol = 0.02;
58 
59 
60 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
61 
62 namespace
63 {
64 
65 // A good ijk index has all components >= 0
66 static inline bool isGoodIndex(const Foam::labelVector& idx)
67 {
68  return (idx.x() >= 0 && idx.y() >= 0 && idx.z() >= 0);
69 }
70 
71 } // End anonymous namespace
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 (
78  const polyMesh& mesh,
79  const PDRblock& pdrBlock
80 )
81 {
82  // Additional copy of i-j-k addressing
83  cellDims = pdrBlock.sizes();
85 
86  const label maxPointId = cmptMax(pdrBlock.sizes())+1;
87 
88  Info<< "Mesh" << nl
89  << " nPoints:" << mesh.nPoints()
90  << " nCells:" << mesh.nCells()
91  << " nFaces:" << mesh.nFaces() << nl;
92 
93  Info<< "PDRblock" << nl
94  << " nPoints:" << pdrBlock.nPoints()
95  << " nCells:" << pdrBlock.nCells()
96  << " nFaces:" << pdrBlock.nFaces() << nl
97  << " min-edge:" << pdrBlock.edgeLimits().min() << nl;
98 
99  Info<< "Classifying ijk indexing... " << nl;
100 
101 
102  // Bin cells into i-j-k locations with the PDRblock::findCell()
103  // method, which combines a bounding box rejection and binary
104  // search in the three directions.
105 
106  cellIndex.resize(mesh.nCells());
107  {
108  const pointField& cc = mesh.cellCentres();
109 
110  for (label celli = 0; celli < mesh.nCells(); ++celli)
111  {
112  cellIndex[celli] = pdrBlock.findCell(cc[celli]);
113  }
114  }
115 
116  // Verify that all i-j-k cells were indeed found
117  {
118  // This could be more efficient - but we want to be picky
119  IjkField<bool> cellFound(pdrBlock.sizes(), false);
120 
121  for (label celli=0; celli < cellIndex.size(); ++celli)
122  {
123  const labelVector& cellIdx = cellIndex[celli];
124 
125  if (isGoodIndex(cellIdx))
126  {
127  cellFound(cellIdx) = true;
128  }
129  }
130 
131  const label firstMiss = cellFound.find(false);
132 
133  if (firstMiss >= 0)
134  {
135  label nMissing = 0;
136  for (label celli = firstMiss; celli < cellFound.size(); ++celli)
137  {
138  if (!cellFound[celli])
139  {
140  ++nMissing;
141  }
142  }
143 
145  << "No ijk location found for "
146  << nMissing << " cells.\nFirst miss at: "
147  << pdrBlock.index(firstMiss)
148  << " indexing" << nl
149  << exit(FatalError);
150  }
151  }
152 
153 
154  // Bin all mesh points into i-j-k locations
155  List<labelVector> pointIndex(mesh.nPoints());
156 
157  for (label pointi = 0; pointi < mesh.nPoints(); ++pointi)
158  {
159  const point& p = mesh.points()[pointi];
160  pointIndex[pointi] = pdrBlock.gridIndex(p, gridPointRelTol);
161  }
162 
163  // Min x,y,z index
164  const labelMinMax invertedLimits(maxPointId, -maxPointId);
165  Vector<labelMinMax> faceLimits;
166 
167  const Vector<direction> faceBits
168  (
172  );
173 
174  faceIndex.resize(mesh.nFaces());
175  faceOrient.resize(mesh.nFaces());
176 
177  for (label facei=0; facei < mesh.nFaces(); ++facei)
178  {
179  labelVector& faceIdx = faceIndex[facei];
180 
181  // Check faces that are associated with i-j-k cells
182  const label own = mesh.faceOwner()[facei];
183  const label nei =
184  (
185  facei < mesh.nInternalFaces()
186  ? mesh.faceNeighbour()[facei]
187  : own
188  );
189 
190  if (!isGoodIndex(cellIndex[own]) && !isGoodIndex(cellIndex[nei]))
191  {
192  // Invalid
193  faceIdx.x() = faceIdx.y() = faceIdx.z() = -1;
194  faceOrient[facei] = vector::X;
195  continue;
196  }
197 
198 
199  faceLimits.x() = faceLimits.y() = faceLimits.z() = invertedLimits;
200 
201  for (const label pointi : mesh.faces()[facei])
202  {
203  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
204  {
205  faceLimits[cmpt].add(pointIndex[pointi][cmpt]);
206  }
207  }
208 
209  direction inPlane(0u);
210 
211  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
212  {
213  const auto& limits = faceLimits[cmpt];
214 
215  if (!limits.valid())
216  {
217  // This should be impossible
219  << "Unexpected search failure for " << facei << " in "
220  << vector::componentNames[cmpt] << "-direction" << nl
221  << exit(FatalError);
222  }
223 
224  if (limits.min() < 0)
225  {
227  << "Face " << facei << " contains non-grid point in "
228  << vector::componentNames[cmpt] << "-direction" << nl
229  << mesh.faces()[facei] << ' '
230  << mesh.faces()[facei].points(mesh.points())
231  << exit(FatalError);
232  }
233  else if (limits.min() == limits.max())
234  {
235  // In plane
236  inPlane |= faceBits[cmpt];
237  }
238  else if (limits.min() + 1 != limits.max())
239  {
241  << "Face " << facei << " not in "
242  << vector::componentNames[cmpt] << "-plane" << nl
243  << exit(FatalError);
244  }
245  }
246 
247  switch (inPlane)
248  {
249  case boundBox::XDIR:
250  faceOrient[facei] = vector::X;
251  break;
252 
253  case boundBox::YDIR:
254  faceOrient[facei] = vector::Y;
255  break;
256 
257  case boundBox::ZDIR:
258  faceOrient[facei] = vector::Z;
259  break;
260 
261  default:
263  << "Face " << facei << " not in an x/y/z plane?" << nl
264  << exit(FatalError);
265  break;
266  }
267 
268  faceIdx.x() = faceLimits.x().min();
269  faceIdx.y() = faceLimits.y().min();
270  faceIdx.z() = faceLimits.z().min();
271  }
272 }
273 
274 
276 (
277  const Time& runTime,
278  const PDRblock& pdrBlock
279 )
280 {
281  #include "createPolyMesh.H"
282  classify(mesh, pdrBlock);
283 }
284 
285 
286 // ************************************************************************* //
Required Variables.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:132
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
labelVector faceDims
The face i-j-k addressing range.
Definition: PDRmeshArrays.H:78
List< labelVector > cellIndex
For each cell, the corresponding i-j-k address.
Definition: PDRmeshArrays.H:83
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
static scalar gridPointRelTol
Relative tolerance when matching grid points. Default = 0.02.
Definition: PDRmeshArrays.H:68
engineTime & runTime
const Cmpt & y() const noexcept
Access to the vector y component.
Definition: Vector.H:140
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
List< labelVector > faceIndex
For each face, the corresponding i-j-k address.
Definition: PDRmeshArrays.H:88
void read(const Time &runTime, const PDRblock &pdrBlock)
Read OpenFOAM mesh and determine i-j-k indices for faces/cells.
List< direction > faceOrient
For each face, the x/y/z orientation.
Definition: PDRmeshArrays.H:93
MinMax< label > labelMinMax
A label min/max range.
Definition: MinMax.H:107
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:109
const Cmpt & x() const noexcept
Access to the vector x component.
Definition: Vector.H:135
4: z-direction. Same as (1 << vector::Z)
Definition: boundBox.H:113
const Cmpt & z() const noexcept
Access to the vector z component.
Definition: Vector.H:145
vector point
Point is a vector.
Definition: point.H:37
PtrList< volScalarField > & Y
1: x-direction. Same as (1 << vector::X)
Definition: boundBox.H:111
messageStream Info
Information stream (stdout output on master, null elsewhere)
void classify(const polyMesh &mesh, const PDRblock &pdrBlock)
Determine i-j-k indices for faces/cells.
Vector< label > labelVector
Vector of labels.
Definition: labelVector.H:47
volScalarField & p
2: y-direction. Same as (1 << vector::Y)
Definition: boundBox.H:112
labelVector cellDims
The cell i-j-k addressing range.
Definition: PDRmeshArrays.H:73