Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / OpenFOAM / meshes / treeBoundBox / treeBoundBoxI.H
blob1f069feb482e0540f82bba0b82126529cdac9354
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "treeBoundBox.H"
27 #include "Random.H"
29 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
31 inline Foam::treeBoundBox::treeBoundBox()
33     boundBox()
37 inline Foam::treeBoundBox::treeBoundBox(const point& min, const point& max)
39     boundBox(min, max)
43 inline Foam::treeBoundBox::treeBoundBox(const boundBox& bb)
45     boundBox(bb)
49 inline Foam::treeBoundBox::treeBoundBox(Istream& is)
51     boundBox(is)
55 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
57 inline Foam::scalar Foam::treeBoundBox::typDim() const
59     return avgDim();
63 inline Foam::point Foam::treeBoundBox::corner(const direction octant) const
65     return point
66     (
67         (octant & RIGHTHALF) ? max().x() : min().x(),
68         (octant & TOPHALF)   ? max().y() : min().y(),
69         (octant & FRONTHALF) ? max().z() : min().z()
70     );
74 // Returns octant in which point resides. Reverse of subBbox.
75 inline Foam::direction Foam::treeBoundBox::subOctant(const point& pt) const
77     return subOctant(midpoint(), pt);
81 // Returns octant in which point resides. Reverse of subBbox.
82 // Precalculated midpoint
83 inline Foam::direction Foam::treeBoundBox::subOctant
85     const point& mid,
86     const point& pt
89     direction octant = 0;
91     if (pt.x() > mid.x())
92     {
93         octant |= treeBoundBox::RIGHTHALF;
94     }
96     if (pt.y() > mid.y())
97     {
98         octant |= treeBoundBox::TOPHALF;
99     }
101     if (pt.z() > mid.z())
102     {
103         octant |= treeBoundBox::FRONTHALF;
104     }
106     return octant;
110 // Returns octant in which point resides. Reverse of subBbox.
111 // Flags point exactly on edge.
112 inline Foam::direction Foam::treeBoundBox::subOctant
114     const point& pt,
115     bool& onEdge
116 ) const
118     return subOctant(midpoint(), pt, onEdge);
122 // Returns octant in which point resides. Reverse of subBbox.
123 // Precalculated midpoint
124 inline Foam::direction Foam::treeBoundBox::subOctant
126     const point& mid,
127     const point& pt,
128     bool& onEdge
131     direction octant = 0;
132     onEdge = false;
134     if (pt.x() > mid.x())
135     {
136         octant |= treeBoundBox::RIGHTHALF;
137     }
138     else if (pt.x() == mid.x())
139     {
140         onEdge = true;
141     }
143     if (pt.y() > mid.y())
144     {
145         octant |= treeBoundBox::TOPHALF;
146     }
147     else if (pt.y() == mid.y())
148     {
149         onEdge = true;
150     }
152     if (pt.z() > mid.z())
153     {
154         octant |= treeBoundBox::FRONTHALF;
155     }
156     else if (pt.z() == mid.z())
157     {
158         onEdge = true;
159     }
161     return octant;
165 // Returns octant in which intersection resides.
166 // Precalculated midpoint. If the point is on the dividing line between
167 // the octants the direction vector determines which octant to use
168 // (i.e. in which octant the point would be if it were moved along dir)
169 inline Foam::direction Foam::treeBoundBox::subOctant
171     const point& mid,
172     const vector& dir,
173     const point& pt,
174     bool& onEdge
177     direction octant = 0;
178     onEdge = false;
180     if (pt.x() > mid.x())
181     {
182         octant |= treeBoundBox::RIGHTHALF;
183     }
184     else if (pt.x() == mid.x())
185     {
186         onEdge = true;
187         if (dir.x() > 0)
188         {
189             octant |= treeBoundBox::RIGHTHALF;
190         }
191     }
193     if (pt.y() > mid.y())
194     {
195         octant |= treeBoundBox::TOPHALF;
196     }
197     else if (pt.y() == mid.y())
198     {
199         onEdge = true;
200         if (dir.y() > 0)
201         {
202             octant |= treeBoundBox::TOPHALF;
203         }
204     }
206     if (pt.z() > mid.z())
207     {
208         octant |= treeBoundBox::FRONTHALF;
209     }
210     else if (pt.z() == mid.z())
211     {
212         onEdge = true;
213         if (dir.z() > 0)
214         {
215             octant |= treeBoundBox::FRONTHALF;
216         }
217     }
219     return octant;
223 // Returns reference to octantOrder which defines the
224 // order to do the search.
225 inline void Foam::treeBoundBox::searchOrder
227     const point& pt,
228     FixedList<direction,8>& octantOrder
229 ) const
231     vector dist = midpoint() - pt;
233     direction octant = 0;
235     if (dist.x() < 0)
236     {
237         octant |= treeBoundBox::RIGHTHALF;
238         dist.x() *= -1;
239     }
241     if (dist.y() < 0)
242     {
243         octant |= treeBoundBox::TOPHALF;
244         dist.y() *= -1;
245     }
247     if (dist.z() < 0)
248     {
249         octant |= treeBoundBox::FRONTHALF;
250         dist.z() *= -1;
251     }
253     direction min = 0;
254     direction mid = 0;
255     direction max = 0;
257     if (dist.x() < dist.y())
258     {
259         if (dist.y() < dist.z())
260         {
261             min = treeBoundBox::RIGHTHALF;
262             mid = treeBoundBox::TOPHALF;
263             max = treeBoundBox::FRONTHALF;
264         }
265         else if (dist.z() < dist.x())
266         {
267             min = treeBoundBox::FRONTHALF;
268             mid = treeBoundBox::RIGHTHALF;
269             max = treeBoundBox::TOPHALF;
270         }
271         else
272         {
273             min = treeBoundBox::RIGHTHALF;
274             mid = treeBoundBox::FRONTHALF;
275             max = treeBoundBox::TOPHALF;
276         }
277     }
278     else
279     {
280         if (dist.z() < dist.y())
281         {
282             min = treeBoundBox::FRONTHALF;
283             mid = treeBoundBox::TOPHALF;
284             max = treeBoundBox::RIGHTHALF;
285         }
286         else if (dist.x() < dist.z())
287         {
288             min = treeBoundBox::TOPHALF;
289             mid = treeBoundBox::RIGHTHALF;
290             max = treeBoundBox::FRONTHALF;
291         }
292         else
293         {
294             min = treeBoundBox::TOPHALF;
295             mid = treeBoundBox::FRONTHALF;
296             max = treeBoundBox::RIGHTHALF;
297         }
298     }
300     // Primary subOctant
301     octantOrder[0] = octant;
302     // subOctants joined to the primary by faces.
303     octantOrder[1] = octant ^ min;
304     octantOrder[2] = octant ^ mid;
305     octantOrder[3] = octant ^ max;
306     // subOctants joined to the primary by edges.
307     octantOrder[4] = octantOrder[1] ^ mid;
308     octantOrder[5] = octantOrder[1] ^ max;
309     octantOrder[6] = octantOrder[2] ^ max;
310     // subOctants joined to the primary by corner.
311     octantOrder[7] = octantOrder[4] ^ max;
315 //- Return slightly wider bounding box
316 inline Foam::treeBoundBox Foam::treeBoundBox::extend
318     Random& rndGen,
319     const scalar s
320 ) const
322     treeBoundBox bb(*this);
324     vector newSpan = bb.span();
326     // Make 3D
327     scalar minSpan = s * Foam::mag(newSpan);
329     for (direction dir = 0; dir < vector::nComponents; dir++)
330     {
331         newSpan[dir] = Foam::max(newSpan[dir], minSpan);
332     }
334     bb.min() -= cmptMultiply(s * rndGen.vector01(), newSpan);
335     bb.max() += cmptMultiply(s * rndGen.vector01(), newSpan);
337     return bb;
341 // ************************************************************************* //