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-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
32{
33 return (scalarList::size() > 1);
34}
35
36
37inline Foam::label Foam::PDRblock::location::nCells() const
38{
39 return (scalarList::size()-1);
40}
41
42
43inline Foam::label Foam::PDRblock::location::nPoints() const
44{
45 return scalarList::size();
46}
47
48
49inline bool Foam::PDRblock::location::contains(const scalar p) const
50{
51 return (scalarList::size() > 1 && first() <= p && p <= last());
52}
53
54
55inline const Foam::scalar& Foam::PDRblock::location::min() const
56{
58}
59
60
61inline const Foam::scalar& Foam::PDRblock::location::max() const
62{
64}
65
66
67inline Foam::scalar Foam::PDRblock::location::centre() const
68{
69 return scalarList::empty() ? 0 : (0.5*first() + 0.5*last());
70}
71
72
73inline Foam::scalar Foam::PDRblock::location::length() const
74{
75 return scalarList::empty() ? 0 : mag(last() - first());
76}
77
78
79inline 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
91inline Foam::scalar Foam::PDRblock::location::width(const label i) const
92{
93 #ifdef FULLDEBUG
94 checkIndex(i);
95 #endif
96
97 return (operator[](i+1) - operator[](i));
98}
99
100
101inline 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 first() - 0.5*(operator[](1) - first());
111 }
112 else if (i > 1 && i == scalarList::size()-1)
113 {
114 // "Halo" centre [nCells] == xN + 1/2 (xN - xN1)
115 return last() + 0.5*(last() - 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
126inline const Foam::scalar&
127Foam::PDRblock::location::clip(const scalar& val) const
128{
129 if (scalarList::size())
130 {
131 if (val < first())
132 {
133 return first();
134 }
135 else if (last() < val)
136 {
137 return last();
138 }
139 }
140
141 return val; // Pass-through
142}
143
144
145// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146
149{
150 return grid_;
151}
152
153
155{
156 return edgeLimits_;
157}
158
159
161{
162 return bounds_;
163}
164
165
166inline Foam::scalar Foam::PDRblock::dx(const label i) const
167{
168 return grid_.x().width(i);
169}
170
171
172inline Foam::scalar Foam::PDRblock::dx(const labelVector& ijk) const
173{
174 return grid_.x().width(ijk.x());
175}
176
177
178inline Foam::scalar Foam::PDRblock::dy(const label j) const
179{
180 return grid_.y().width(j);
181}
182
183
184inline Foam::scalar Foam::PDRblock::dy(const labelVector& ijk) const
185{
186 return grid_.y().width(ijk.y());
187}
188
189
190inline Foam::scalar Foam::PDRblock::dz(const label k) const
191{
192 return grid_.z().width(k);
193}
194
195
196inline Foam::scalar Foam::PDRblock::dz(const labelVector& ijk) const
197{
198 return grid_.z().width(ijk.z());
199}
200
201
203(
204 const label i,
205 const label j,
206 const label k
207) const
208{
209 return vector(dx(i), dy(j), dz(k));
210}
211
212
214{
215 return vector(dx(ijk), dy(ijk), dz(ijk));
216}
217
218
220(
221 const label i,
222 const label j,
223 const label k
224) const
225{
226 return point(grid_.x()[i], grid_.y()[j], grid_.z()[k]);
227}
228
229
231{
232 return
233 point
234 (
235 grid_.x()[ijk.x()],
236 grid_.y()[ijk.y()],
237 grid_.z()[ijk.z()]
238 );
239}
240
241
243(
244 const label i,
245 const label j,
246 const label k
247) const
248{
249 return point(grid_.x().C(i), grid_.y().C(j), grid_.z().C(k));
250}
251
252
254{
255 return
256 point
257 (
258 grid_.x().C(ijk.x()),
259 grid_.y().C(ijk.y()),
260 grid_.z().C(ijk.z())
261 );
262}
263
264
265inline Foam::scalar Foam::PDRblock::V
266(
267 const label i,
268 const label j,
269 const label k
270) const
271{
272 return dx(i)*dy(j)*dz(k);
273}
274
275
276inline Foam::scalar Foam::PDRblock::V(const labelVector& ijk) const
277{
278 return dx(ijk.x())*dy(ijk.y())*dz(ijk.z());
279}
280
281
282inline Foam::scalar Foam::PDRblock::width
283(
284 const label i,
285 const label j,
286 const label k
287) const
288{
289 return Foam::cbrt(V(i, j, k));
290}
291
292
293inline Foam::scalar Foam::PDRblock::width(const labelVector& ijk) const
294{
295 return Foam::cbrt(V(ijk));
296}
297
298
299// ************************************************************************* //
label k
Graphite solid properties.
Definition: C.H:53
scalar length() const
The difference between min/max values, zero for an empty list.
Definition: PDRblockI.H:73
bool valid() const
The location list is valid if it contains 2 or more points.
Definition: PDRblockI.H:31
void checkIndex(const label i) const
Check that element index is within valid range.
Definition: PDRblockI.H:79
scalar centre() const
Mid-point location, zero for an empty list.
Definition: PDRblockI.H:67
const scalar & min() const
The first() value is considered the min value.
Definition: PDRblockI.H:55
label nCells() const
The number of cells in this direction.
Definition: PDRblockI.H:37
label nPoints() const
The number of points in this direction.
Definition: PDRblockI.H:43
const scalar & max() const
The last() value is considered the max value.
Definition: PDRblockI.H:61
bool contains(const scalar p) const
True if the location is within the range.
Definition: PDRblockI.H:49
const scalar & clip(const scalar &val) const
Definition: PDRblockI.H:127
scalar V(const label i, const label j, const label k) const
Cell volume at i,j,k position.
Definition: PDRblockI.H:266
const boundBox & bounds() const
The mesh bounding box.
Definition: PDRblockI.H:160
scalar dx(const label i) const
Cell size in x-direction at i position.
Definition: PDRblockI.H:166
const scalarMinMax & edgeLimits() const
The min/max edge length.
Definition: PDRblockI.H:154
scalar dy(const label j) const
Cell size in y-direction at j position.
Definition: PDRblockI.H:178
scalar dz(const label k) const
Cell size in z-direction at k position.
Definition: PDRblockI.H:190
const Vector< location > & grid() const
The grid point locations in the i,j,k (x,y,z) directions.
Definition: PDRblockI.H:148
vector span
The obstacle dimensions (for boxes)
Definition: PDRobstacle.H:126
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
label width() const noexcept
Return current width of mask and padded.
Definition: ensightCase.H:428
void checkIndex(const label i, const label j, const label k, const bool allowExtra=false) const
Check indices are within ni,nj,nk range.
label nCells() const
The number of mesh cells (nx*ny*nz) in the i-j-k mesh.
Definition: ijkMeshI.H:68
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
scalar dx() const
Return the fixed interval.
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
vector point
Point is a vector.
Definition: point.H:43
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
Vector< scalar > vector
Definition: vector.H:61
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53