PDRutilsInternal.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-------------------------------------------------------------------------------
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
26Namespace
27 Foam::PDRutils
28
29Description
30 Utilities for PDR (eg, for setFields). Internal usage only.
31
32 The C lineage of the original code is still evident in the use of
33 pointers instead of references.
34 This will be addressed in later versions of the code (2019-12).
35
36SourceFiles
37 PDRUtils.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef PDRutilsInternal_H
42#define PDRutilsInternal_H
43
44#include "PDRutils.H"
45#include "PDRarrays.H"
46#include "PDRblock.H"
47#include "symmTensor2D.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53namespace PDRutils
54{
55
56//- Determine 1-D overlap locations for a geometric entity
57//
58// \param[in] xmin - min position of the geometric entity
59// \param[in] xmax - max position of the geometric entity
60// \param[in] grid - grid point information
61// \param[out] olap - Fraction of cell-width with overlap
62// 0 for no overlap, 1 for full overlap.
63// \param[out] cmin - first cell index (inclusive) with overlap,
64// values in the range \c [0,nCells]
65// \param[out] cmax - last cell index (inclusive) with overlap,
66// values in the range \c [0,nCells]
67// \param[out] cfmin - first cell index (inclusive) with staggered face,
68// values in the range \c [0,nCells]
69// \param[out] cfmax - last cell index (inclusive) with staggered face,
70// values in the range \c [0,nCells]
72(
73 scalar xmin,
74 scalar xmax,
75 const PDRblock::location& grid,
76 List<scalar>& olap,
77 int *cmin, int *cmax,
78 int *cfmin, int *cfmax
79);
80
81
82//- Combine two 1D overlaps.
83// Multiplying the two 1-d overlaps yields the proportion of each (2D) cell
84// that is covered.
85//
86// \note We go one over the relevant min/max limits since these
87// values might be used.
88// The 1D arrays will have bee initially zeroed throughout.
90(
91 const UList<scalar>& a_olap, label amin, label amax,
92 const UList<scalar>& b_olap, label bmin, label bmax,
94);
95
96
97//- Calculate the proportion of each (two-dimensional) grid cell
98//- overlapped by the circle or angled rectangle.
99//
100// Coordinates are labelled a and b.
101//
102// \param[in] ac, bc coordinates of centre of circle or rectangle
103// \param[in] dia diameter of circle (zero for rectangle)
104// \param[in] theta, wa, wb parameters for rectangle
105// \param[in] amin, amax first and last cells in a-grid overlapped by object
106// \param[in] agrid locations of grid lines of a-grid
107// \param[in] amin, amax first and last cells in b-grid overlapped by object
108// \param[in] bgrid locations of grid lines of b-grid
109//
110// \param[out] abolap
111// 2-D array of (proportionate) area blockage by grid cell
112// \param[out] a_lblock
113// 2-D array of (proportionate) blockage to a-direction flow
114// (This will be area blockage when extruded in the third
115// coordinate).
116//
117// \param[out] a_count
118// 2-D array The contribution of this object to the count of
119// obstacles blocking a-direction flow. This is only non-zero if the
120// object is inside the lateral boundaries of the cell. It is large
121// negative if the cell is totally blocked in this direction.
122//
123//
124// \param[out] c_drag
125//
126// 2-D array of tensor that will give tensor drag in each cell (when
127// multiplied Cd, cylinder length, and 0.5 rho*U^2) Dimension: L.
128//
129// \note this routine does not zero array elements outside the amin
130// to amax, bmin to bmax area.
132(
133 scalar ac, scalar bc, scalar dia,
134 scalar theta, scalar wa, scalar wb,
135 const PDRblock::location& agrid, label amin, label amax,
136 const PDRblock::location& bgrid, label bmin, label bmax,
137 SquareMatrix<scalar>& ab_olap,
138 SquareMatrix<scalar>& ab_perim,
139 SquareMatrix<scalar>& a_lblock,
140 SquareMatrix<scalar>& ac_lblock,
141 SquareMatrix<scalar>& c_count,
143 SquareMatrix<scalar>& b_lblock,
144 SquareMatrix<scalar>& bc_lblock
145);
146
147
148//- Area of intersection between circle and rectangle.
149//
150// Calculates the area of intersection between the circle, centre (xc, yc), radius rad,
151// and the rectangle with sides at x = x1 & x2, and y = y1 and y2.
152//
153// The return value is the fraction of the rectangle's area covered by the circle.
155(
156 double xc,
157 double yc,
158 double rad,
159 double x1, double x2,
160 double y1, double y2,
161 scalar* perim_p,
162 scalar* x_proj_edge_p,
163 scalar* y_proj_edge_p,
164 scalar* x_overlap_p,
165 scalar* y_overlap_p
166);
167
168
169//- The area overlap in the plane of a diagonal block and a cell.
170//
171// Calculates the overlap, in the plane of a diagonal block and a cell,
172// plus blockage and drag parameters.
173// Note that x and y herein may be any two of the three coordinates - would have been better not to label them x and y.
174//
175// On entry:
176// xc, yc Coordinates of axis of d.b.
177// theta, wa, wb Angle and widths
178//
179// The returned parameters will be multiplied by the length of the obstacle's intersection with
180// the third dimension of the 3-D cell to give this obstacle's contribution to the count, drag
181// and area blockages.
182// The return value is the area of intersection, which will multiply to volume blockage.
183//
185(
186 double xc, double yc, double theta,
187 double wa, double wb,
188 double x1, double x2,
189 double y1, double y2,
190 scalar* count_p,
191 symmTensor2D& vdrag, scalar* perim_p,
192 scalar* x_lblk, scalar* y_lblk,
193 scalar* x_centre_p, scalar* y_centre_p
194);
195
196
197/* Calculates the blockage to x-direction flow presented by the specified circle on
198 the specified rectangle.
199 Becomes the area blockage when extruded to in the third dimension.
200 In other words, it is the projection on the y axis of the intersection between the
201 circle and the rectangle.
202 Returns fraction blocked
203 Note that x and y in this routine may in fact be any two of the three dimensions.
204 */
206(
207 double xc, double yc, double rad,
208 double x1, double x2,
209 double y1, double y2,
210 scalar* count_p, scalar* drag_p, scalar* centre_p
211);
212
213
214} // End namespace PDRutils
215} // End namespace Foam
216
217// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218
219#endif
220
221// ************************************************************************* //
radiation::radiationModel & rad
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
Grid locations in an axis direction.
Definition: PDRblock.H:181
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
Definition: SquareMatrix.H:66
A templated (2 x 2) symmetric tensor of objects of <T>, effectively containing 3 elements,...
Definition: SymmTensor2D.H:61
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void two_d_overlap(const UList< scalar > &a_olap, label amin, label amax, const UList< scalar > &b_olap, label bmin, label bmax, SquareMatrix< scalar > &ab_olap)
Combine two 1D overlaps.
void one_d_overlap(scalar xmin, scalar xmax, const PDRblock::location &grid, List< scalar > &olap, int *cmin, int *cmax, int *cfmin, int *cfmax)
Determine 1-D overlap locations for a geometric entity.
void circle_overlap(scalar ac, scalar bc, scalar dia, scalar theta, scalar wa, scalar wb, const PDRblock::location &agrid, label amin, label amax, const PDRblock::location &bgrid, label bmin, label bmax, SquareMatrix< scalar > &ab_olap, SquareMatrix< scalar > &ab_perim, SquareMatrix< scalar > &a_lblock, SquareMatrix< scalar > &ac_lblock, SquareMatrix< scalar > &c_count, SquareMatrix< symmTensor2D > &c_drag, SquareMatrix< scalar > &b_lblock, SquareMatrix< scalar > &bc_lblock)
double l_blockage(double xc, double yc, double rad, double x1, double x2, double y1, double y2, scalar *count_p, scalar *drag_p, scalar *centre_p)
double inters_db(double xc, double yc, double theta, double wa, double wb, double x1, double x2, double y1, double y2, scalar *count_p, symmTensor2D &vdrag, scalar *perim_p, scalar *x_lblk, scalar *y_lblk, scalar *x_centre_p, scalar *y_centre_p)
The area overlap in the plane of a diagonal block and a cell.
double inters_cy(double xc, double yc, double rad, double x1, double x2, double y1, double y2, scalar *perim_p, scalar *x_proj_edge_p, scalar *y_proj_edge_p, scalar *x_overlap_p, scalar *y_overlap_p)
Area of intersection between circle and rectangle.
Namespace for OpenFOAM.
SymmTensor2D< scalar > symmTensor2D
SymmTensor2D of scalars, i.e. SymmTensor2D<scalar>.
Definition: symmTensor2D.H:66
dimensionedScalar y1(const dimensionedScalar &ds)