PDRblockI.H
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) 2019 OpenCFD Ltd.
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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 :
32  ijkMesh(),
33  verbose_(false),
34  grid_(),
35  bounds_(),
36  patches_(),
37  minEdgeLen_(Zero)
38 {}
39 
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
44 {
45  return (scalarList::size() > 1);
46 }
47 
48 
50 {
51  return (scalarList::size()-1);
52 }
53 
54 
56 {
57  return scalarList::size();
58 }
59 
60 
61 inline bool Foam::PDRblock::location::contains(const scalar p) const
62 {
63  return (scalarList::size() > 1 && first() <= p && p <= last());
64 }
65 
66 
67 inline const Foam::scalar& Foam::PDRblock::location::min() const
68 {
69  return scalarList::empty() ? pTraits<scalar>::rootMax : first();
70 }
71 
72 
73 inline const Foam::scalar& Foam::PDRblock::location::max() const
74 {
75  return scalarList::empty() ? pTraits<scalar>::rootMin : last();
76 }
77 
78 
79 inline Foam::scalar Foam::PDRblock::location::centre() const
80 {
81  return scalarList::empty() ? 0 : (0.5*first() + 0.5*last());
82 }
83 
84 
85 inline void Foam::PDRblock::location::checkIndex(const label i) const
86 {
87  if (i < 0 || i >= nCells())
88  {
90  << "The index " << i
91  << " is out of range [0," << nCells() << ']' << nl
92  << abort(FatalError);
93  }
94 }
95 
96 
97 inline Foam::scalar Foam::PDRblock::location::width(const label i) const
98 {
99  #ifdef FULLDEBUG
100  checkIndex(i);
101  #endif
102 
103  return (operator[](i+1) - operator[](i));
104 }
105 
106 
107 inline Foam::scalar Foam::PDRblock::location::C(const label i) const
108 {
109  if (i == -1)
110  {
111  #ifdef FULLDEBUG
112  checkIndex(0);
113  #endif
114 
115  // "Halo" centre [-1] == x0 - 1/2 (x1 - x0)
116  return first() - 0.5*(operator[](1) - first());
117  }
118  else if (i > 1 && i == scalarList::size()-1)
119  {
120  // "Halo" centre [nCells] == xN + 1/2 (xN - xN1)
121  return last() + 0.5*(last() - operator[](scalarList::size()-2));
122  }
123 
124  #ifdef FULLDEBUG
125  checkIndex(i);
126  #endif
127 
128  return 0.5*(operator[](i+1) + operator[](i));
129 }
130 
131 
132 inline const Foam::scalar&
134 {
135  if (scalarList::size())
136  {
137  if (val < first())
138  {
139  return first();
140  }
141  else if (last() < val)
142  {
143  return last();
144  }
145  }
146 
147  return val; // Pass-through
148 }
149 
150 
153 {
154  return grid_;
155 }
156 
157 
158 inline const Foam::scalar& Foam::PDRblock::minEdgeLen() const
159 {
160  return minEdgeLen_;
161 }
162 
163 
165 {
166  return bounds_;
167 }
168 
169 
170 inline Foam::scalar Foam::PDRblock::dx(const label i) const
171 {
172  return grid_.x().width(i);
173 }
174 
175 
176 inline Foam::scalar Foam::PDRblock::dx(const labelVector& ijk) const
177 {
178  return grid_.x().width(ijk.x());
179 }
180 
181 
182 inline Foam::scalar Foam::PDRblock::dy(const label j) const
183 {
184  return grid_.y().width(j);
185 }
186 
187 
188 inline Foam::scalar Foam::PDRblock::dy(const labelVector& ijk) const
189 {
190  return grid_.y().width(ijk.y());
191 }
192 
193 
194 inline Foam::scalar Foam::PDRblock::dz(const label k) const
195 {
196  return grid_.z().width(k);
197 }
198 
199 
200 inline Foam::scalar Foam::PDRblock::dz(const labelVector& ijk) const
201 {
202  return grid_.z().width(ijk.z());
203 }
204 
205 
207 (
208  const label i,
209  const label j,
210  const label k
211 ) const
212 {
213  return vector(dx(i), dy(j), dz(k));
214 }
215 
216 
218 {
219  return vector(dx(ijk), dy(ijk), dz(ijk));
220 }
221 
222 
224 (
225  const label i,
226  const label j,
227  const label k
228 ) const
229 {
230  return point(grid_.x()[i], grid_.y()[j], grid_.z()[k]);
231 }
232 
233 
235 {
236  return
237  point
238  (
239  grid_.x()[ijk.x()],
240  grid_.y()[ijk.y()],
241  grid_.z()[ijk.z()]
242  );
243 }
244 
245 
247 (
248  const label i,
249  const label j,
250  const label k
251 ) const
252 {
253  return point(grid_.x().C(i), grid_.y().C(j), grid_.z().C(k));
254 }
255 
256 
257 inline Foam::point Foam::PDRblock::C(const labelVector& ijk) const
258 {
259  return
260  point
261  (
262  grid_.x().C(ijk.x()),
263  grid_.y().C(ijk.y()),
264  grid_.z().C(ijk.z())
265  );
266 }
267 
268 
269 inline Foam::scalar Foam::PDRblock::V
270 (
271  const label i,
272  const label j,
273  const label k
274 ) const
275 {
276  return dx(i)*dy(j)*dz(k);
277 }
278 
279 
280 inline Foam::scalar Foam::PDRblock::V(const labelVector& ijk) const
281 {
282  return dx(ijk.x())*dy(ijk.y())*dz(ijk.z());
283 }
284 
285 
286 inline Foam::scalar Foam::PDRblock::width
287 (
288  const label i,
289  const label j,
290  const label k
291 ) const
292 {
293  return Foam::cbrt(V(i, j, k));
294 }
295 
296 
297 inline Foam::scalar Foam::PDRblock::width(const labelVector& ijk) const
298 {
299  return Foam::cbrt(V(ijk));
300 }
301 
302 
303 // ************************************************************************* //
Foam::PDRblock::width
scalar width(const label i, const label j, const label k) const
Characteristic cell size at i,j,k position.
Definition: PDRblockI.H:287
Foam::PDRblock::dz
scalar dz(const label k) const
Cell size in z-direction at k position.
Definition: PDRblockI.H:194
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:78
Foam::val
label ListType::const_reference val
Definition: ListOps.H:407
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::PDRblock::location::valid
bool valid() const
The locations are valid if they contain 2 or more points.
Definition: PDRblockI.H:43
Foam::PDRblock::dy
scalar dy(const label j) const
Cell size in y-direction at j position.
Definition: PDRblockI.H:182
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::PDRblock::location::contains
bool contains(const scalar p) const
True if the location is within the range.
Definition: PDRblockI.H:61
Foam::PDRblock::minEdgeLen
const scalar & minEdgeLen() const
The min edge length.
Definition: PDRblockI.H:158
Foam::PDRblock::location::checkIndex
void checkIndex(const label i) const
Check that element index is within valid range.
Definition: PDRblockI.H:85
Foam::PDRblock::location::min
const scalar & min() const
The first() value is considered the min value.
Definition: PDRblockI.H:67
Foam::PDRblock::bounds
const boundBox & bounds() const
The mesh bounding box.
Definition: PDRblockI.H:164
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:90
Foam::PDRblock::location::clip
const scalar & clip(const scalar &val) const
Definition: PDRblockI.H:133
Foam::PDRblock::dx
scalar dx(const label i) const
Cell size in x-direction at i position.
Definition: PDRblockI.H:170
Foam::label
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:62
Foam::ijkMesh
A simple i-j-k (row-major order) to linear addressing for a rectilinear mesh. Since the underlying me...
Definition: ijkMesh.H:53
Foam::PDRblock::C
point C(const label i, const label j, const label k) const
Cell centre at i,j,k position.
Definition: PDRblockI.H:247
Foam::PDRblock::location::centre
scalar centre() const
Mid-point location, zero for an empty list.
Definition: PDRblockI.H:79
Foam::FatalError
error FatalError
Foam::PDRblock::V
scalar V(const label i, const label j, const label k) const
Cell volume at i,j,k position.
Definition: PDRblockI.H:270
Foam::PDRblock::grid
const Vector< location > & grid() const
The grid point locations in the i,j,k (x,y,z) directions.
Definition: PDRblockI.H:152
Foam::PDRblock::location::nPoints
label nPoints() const
The number of points in this direction.
Definition: PDRblockI.H:55
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::PDRblock::location::width
scalar width(const label i) const
Cell size at element position.
Definition: PDRblockI.H:97
Foam::PDRblock::span
vector span(const label i, const label j, const label k) const
Cell dimensions at i,j,k position.
Definition: PDRblockI.H:207
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:84
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::Vector< Foam::PDRblock::location >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::PDRblock::location::nCells
label nCells() const
The number of cells in this direction.
Definition: PDRblockI.H:49
Foam::PDRblock::location::max
const scalar & max() const
The last() value is considered the max value.
Definition: PDRblockI.H:73
Foam::PDRblock::PDRblock
PDRblock()
Construct zero-size.
Definition: PDRblockI.H:30
Foam::PDRblock::location::C
scalar C(const label i) const
Cell centre at element position.
Definition: PDRblockI.H:107