1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2010-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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 \*---------------------------------------------------------------------------*/
27 #include "pointField.H"
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 inline Foam::boundBox::boundBox()
39 inline Foam::boundBox::boundBox(const point& min, const point& max)
46 inline Foam::boundBox::boundBox(Istream& is)
48 operator>>(is, *this);
52 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 inline const Foam::point& Foam::boundBox::min() const
60 inline const Foam::point& Foam::boundBox::max() const
66 inline Foam::point& Foam::boundBox::min()
72 inline Foam::point& Foam::boundBox::max()
78 inline Foam::point Foam::boundBox::midpoint() const
80 return 0.5 * (max_ + min_);
84 inline Foam::vector Foam::boundBox::span() const
90 inline Foam::scalar Foam::boundBox::mag() const
92 return ::Foam::mag(max_ - min_);
96 inline Foam::scalar Foam::boundBox::volume() const
98 return cmptProduct(span());
102 inline Foam::scalar Foam::boundBox::minDim() const
104 return cmptMin(span());
108 inline Foam::scalar Foam::boundBox::maxDim() const
110 return cmptMax(span());
114 inline Foam::scalar Foam::boundBox::avgDim() const
116 return cmptAv(span());
120 inline bool Foam::boundBox::overlaps(const boundBox& bb) const
124 bb.max_.x() >= min_.x() && bb.min_.x() <= max_.x()
125 && bb.max_.y() >= min_.y() && bb.min_.y() <= max_.y()
126 && bb.max_.z() >= min_.z() && bb.min_.z() <= max_.z()
131 inline bool Foam::boundBox::contains(const point& pt) const
135 pt.x() >= min_.x() && pt.x() <= max_.x()
136 && pt.y() >= min_.y() && pt.y() <= max_.y()
137 && pt.z() >= min_.z() && pt.z() <= max_.z()
142 // this.bb fully contains bb
143 inline bool Foam::boundBox::contains(const boundBox& bb) const
145 return contains(bb.min()) && contains(bb.max());
149 inline bool Foam::boundBox::containsInside(const point& pt) const
153 pt.x() > min_.x() && pt.x() < max_.x()
154 && pt.y() > min_.y() && pt.y() < max_.y()
155 && pt.z() > min_.z() && pt.z() < max_.z()
160 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
162 inline bool Foam::operator==(const boundBox& a, const boundBox& b)
164 return (a.min_ == b.min_) && (a.max_ == b.max_);
168 inline bool Foam::operator!=(const boundBox& a, const boundBox& b)
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //