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-2023 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 // Location
30 
31 inline bool Foam::PDRblock::location::good() const noexcept
32 {
33  return (scalarList::size() > 1);
34 }
35 
36 
37 inline Foam::label Foam::PDRblock::location::nCells() const noexcept
38 {
39  return (scalarList::size()-1);
40 }
41 
42 
43 inline Foam::label Foam::PDRblock::location::nPoints() const noexcept
44 {
45  return scalarList::size();
46 }
47 
48 
49 inline bool Foam::PDRblock::location::contains(const scalar p) const
50 {
51  return (scalarList::size() > 1 && front() <= p && p <= back());
52 }
53 
54 
55 inline const Foam::scalar& Foam::PDRblock::location::min() const
56 {
57  return scalarList::empty() ? pTraits<scalar>::rootMax : front();
58 }
59 
60 
61 inline const Foam::scalar& Foam::PDRblock::location::max() const
62 {
63  return scalarList::empty() ? pTraits<scalar>::rootMin : back();
64 }
65 
66 
67 inline Foam::scalar Foam::PDRblock::location::centre() const
68 {
69  return scalarList::empty() ? 0 : (0.5*front() + 0.5*back());
70 }
71 
72 
73 inline Foam::scalar Foam::PDRblock::location::length() const
74 {
75  return scalarList::empty() ? 0 : mag(back() - front());
76 }
77 
78 
79 inline void Foam::PDRblock::location::checkIndex(const label i) const
80 {
81  if (i < 0 || i >= nCells())
82  {
84  << "The index " << i
85  << " is out of range [0," << nCells() << ']' << nl
86  << abort(FatalError);
87  }
88 }
89 
90 
91 inline Foam::scalar Foam::PDRblock::location::width(const label i) const
92 {
93  #ifdef FULLDEBUG
95  #endif
96 
97  return (operator[](i+1) - operator[](i));
98 }
99 
100 
101 inline Foam::scalar Foam::PDRblock::location::C(const label i) const
102 {
103  if (i == -1)
104  {
105  #ifdef FULLDEBUG
106  checkIndex(0);
107  #endif
108 
109  // "Halo" centre [-1] == x0 - 1/2 (x1 - x0)
110  return front() - 0.5*(operator[](1) - front());
111  }
112  else if (i > 1 && i == scalarList::size()-1)
113  {
114  // "Halo" centre [nCells] == xN + 1/2 (xN - xN1)
115  return back() + 0.5*(back() - operator[](scalarList::size()-2));
116  }
117 
118  #ifdef FULLDEBUG
119  checkIndex(i);
120  #endif
121 
122  return 0.5*(operator[](i+1) + operator[](i));
123 }
124 
125 
126 inline const Foam::scalar&
127 Foam::PDRblock::location::clip(const scalar& val) const
128 {
129  if (scalarList::size())
130  {
131  if (val < front())
132  {
133  return front();
134  }
135  else if (back() < val)
136  {
137  return back();
138  }
139  }
141  return val; // Pass-through
142 }
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 inline Foam::scalar Foam::PDRblock::dx(const label i) const
148 {
149  return grid_.x().width(i);
150 }
151 
153 inline Foam::scalar Foam::PDRblock::dx(const labelVector& ijk) const
154 {
155  return grid_.x().width(ijk.x());
156 }
157 
159 inline Foam::scalar Foam::PDRblock::dy(const label j) const
160 {
161  return grid_.y().width(j);
162 }
163 
165 inline Foam::scalar Foam::PDRblock::dy(const labelVector& ijk) const
166 {
167  return grid_.y().width(ijk.y());
168 }
169 
171 inline Foam::scalar Foam::PDRblock::dz(const label k) const
172 {
173  return grid_.z().width(k);
174 }
175 
176 
177 inline Foam::scalar Foam::PDRblock::dz(const labelVector& ijk) const
178 {
179  return grid_.z().width(ijk.z());
180 }
181 
182 
184 (
185  const label i,
186  const label j,
187  const label k
188 ) const
189 {
190  return vector(dx(i), dy(j), dz(k));
191 }
192 
193 
195 {
196  return vector(dx(ijk), dy(ijk), dz(ijk));
197 }
198 
199 
201 (
202  const label i,
203  const label j,
204  const label k
205 ) const
206 {
207  return point(grid_.x()[i], grid_.y()[j], grid_.z()[k]);
208 }
209 
210 
211 inline Foam::point Foam::PDRblock::grid(const labelVector& ijk) const
212 {
213  return
214  point
215  (
216  grid_.x()[ijk.x()],
217  grid_.y()[ijk.y()],
218  grid_.z()[ijk.z()]
219  );
220 }
221 
222 
224 (
225  const label i,
226  const label j,
227  const label k
228 ) const
229 {
230  return point(grid_.x().C(i), grid_.y().C(j), grid_.z().C(k));
231 }
232 
233 
234 inline Foam::point Foam::PDRblock::C(const labelVector& ijk) const
235 {
236  return
237  point
238  (
239  grid_.x().C(ijk.x()),
240  grid_.y().C(ijk.y()),
241  grid_.z().C(ijk.z())
242  );
243 }
244 
245 
246 inline Foam::scalar Foam::PDRblock::V
247 (
248  const label i,
249  const label j,
250  const label k
251 ) const
252 {
253  return dx(i)*dy(j)*dz(k);
254 }
255 
256 
257 inline Foam::scalar Foam::PDRblock::V(const labelVector& ijk) const
258 {
259  return dx(ijk.x())*dy(ijk.y())*dz(ijk.z());
260 }
261 
262 
263 inline Foam::scalar Foam::PDRblock::width
264 (
265  const label i,
266  const label j,
267  const label k
268 ) const
269 {
270  return Foam::cbrt(V(i, j, k));
271 }
272 
273 
274 inline Foam::scalar Foam::PDRblock::width(const labelVector& ijk) const
275 {
276  return Foam::cbrt(V(ijk));
277 }
278 
279 
280 // ************************************************************************* //
scalar centre() const
Mid-point location, zero for an empty list.
Definition: PDRblockI.H:60
const scalar & max() const
The last() value is considered the max value.
Definition: PDRblockI.H:54
point C(const label i, const label j, const label k) const
Cell centre at i,j,k position.
Definition: PDRblockI.H:217
const scalar & clip(const scalar &val) const
If out of range, return the respective min/max limits, otherwise return the value itself...
Definition: PDRblockI.H:120
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
scalar V(const label i, const label j, const label k) const
Cell volume at i,j,k position.
Definition: PDRblockI.H:240
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void checkIndex(const label i) const
Check that element index is within valid range.
Definition: PDRblockI.H:72
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
label nCells() const
The number of mesh cells (nx*ny*nz) in the i-j-k mesh.
Definition: ijkMeshI.H:61
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
const scalar & min() const
The first() value is considered the min value.
Definition: PDRblockI.H:48
vector span(const label i, const label j, const label k) const
Cell dimensions at i,j,k position.
Definition: PDRblockI.H:177
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
const Cmpt & y() const noexcept
Access to the vector y component.
Definition: Vector.H:140
label k
Boltzmann constant.
label nCells() const noexcept
The number of cells in this direction.
Definition: PDRblockI.H:30
bool contains(const scalar p) const
True if the location is within the range.
Definition: PDRblockI.H:42
bool good() const noexcept
The location list is valid if it contains 2 or more points.
Definition: PDRblockI.H:24
scalar dz(const label k) const
Cell size in z-direction at k position.
Definition: PDRblockI.H:164
dimensionedScalar cbrt(const dimensionedScalar &ds)
Vector< scalar > vector
Definition: vector.H:57
scalar width(const label i, const label j, const label k) const
Characteristic cell size at i,j,k position.
Definition: PDRblockI.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void checkIndex(const label i, const label j, const label k, const bool allowExtra=false) const
Check indices are within ni,nj,nk range.
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
const Cmpt & x() const noexcept
Access to the vector x component.
Definition: Vector.H:135
scalar length() const
The difference between min/max values, zero for an empty list.
Definition: PDRblockI.H:66
scalar dy(const label j) const
Cell size in y-direction at j position.
Definition: PDRblockI.H:152
scalar width(const label i) const
Cell size at element position.
Definition: PDRblockI.H:84
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
volScalarField & p
scalar C(const label i) const
Cell centre at element position.
Definition: PDRblockI.H:94
scalar dx(const label i) const
Cell size in x-direction at i position.
Definition: PDRblockI.H:140
label nPoints() const noexcept
The number of points in this direction.
Definition: PDRblockI.H:36
const Vector< location > & grid() const noexcept
The grid point locations in the i,j,k (x,y,z) directions.
Definition: PDRblock.H:711