FreeFOAM The Cross-Platform CFD Toolkit
PointHit_.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::PointHit
26 
27 Description
28  This class describes the interaction of a face and a point. It
29  carries the info of a successful hit and (if successful), returns
30  the interaction point.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef PointHit__H
35 #define PointHit__H
36 
37 #include <OpenFOAM/bool.H>
38 #include <OpenFOAM/token.H>
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 // Forward declaration of classes
46 
47 class Ostream;
48 
49 
50 // Forward declaration of friend functions and operators
51 
52 template<class Point> class PointHit;
53 
54 template<class Point>
55 inline Ostream& operator<<(Ostream&, const PointHit<Point>&);
56 
57 
58 /*---------------------------------------------------------------------------*\
59  Class PointHit Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 template<class Point>
63 class PointHit
64 {
65  // Private data
66 
67  //- Hit success
68  bool hit_;
69 
70  //- Point of hit; for miss holds best estimate outside the object
71  Point hitPoint_;
72 
73  //- Distance to hit point
74  scalar distance_;
75 
76  //- Eligible miss
77  bool eligibleMiss_;
78 
79 
80 public:
81 
82  // Constructors
83 
84  //- Construct from components
85  PointHit
86  (
87  const bool hit,
88  const Point& p,
89  const scalar dist,
90  const bool eligibleMiss
91  )
92  :
93  hit_(hit),
94  hitPoint_(p),
95  distance_(dist),
96  eligibleMiss_(eligibleMiss)
97  {}
98 
99  //- Construct from point. Hit and distance set later
100  PointHit(const Point& p)
101  :
102  hit_(false),
103  hitPoint_(p),
104  distance_(GREAT),
105  eligibleMiss_(false)
106  {}
107 
108 
109  // Member Functions
110 
111  //- Is there a hit
112  bool hit() const
113  {
114  return hit_;
115  }
116 
117  //- Return hit point
118  const Point& hitPoint() const
119  {
120  if (!hit_)
121  {
122  FatalErrorIn("const Point& PointHit::hitPoint() const")
123  << "requested a hit point for a miss"
124  << abort(FatalError);
125  }
126 
127  return hitPoint_;
128  }
129 
130  //- Return distance to hit
131  scalar distance() const
132  {
133  return distance_;
134  }
135 
136  //- Return miss point
137  const Point& missPoint() const
138  {
139  if (hit_)
140  {
141  FatalErrorIn("const Point& PointHit::missPoint() const")
142  << "requested a miss point for a hit"
143  << abort(FatalError);
144  }
145 
146  return hitPoint_;
147  }
148 
149  //- Return point with no checking
150  const Point& rawPoint() const
151  {
152  return hitPoint_;
153  }
154 
155  //- Is this an eligible miss
156  bool eligibleMiss() const
157  {
158  return eligibleMiss_;
159  }
160 
161  void setHit()
162  {
163  hit_ = true;
164  eligibleMiss_ = false;
165  }
166 
167  void setMiss(const bool eligible)
168  {
169  hit_ = false;
170  eligibleMiss_ = eligible;
171  }
172 
173  void setPoint(const Point& p)
174  {
175  hitPoint_ = p;
176  }
177 
178  void setDistance(const scalar d)
179  {
180  distance_ = d;
181  }
182 
183 
184  // Ostream operator
185 
186  friend Ostream& operator<< <Point>
187  (
188  Ostream& os,
189  const PointHit<Point>& b
190  );
191 };
192 
193 
194 template<class Point>
195 inline Ostream& operator<<(Ostream& os, const PointHit<Point>& b)
196 {
197  os << b.hit() << token::SPACE
198  << b.rawPoint() << token::SPACE
199  << b.distance() << token::SPACE
200  << b.eligibleMiss();
201 
202  return os;
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 } // End namespace Foam
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 #endif
213 
214 // ************************ vim: set sw=4 sts=4 et: ************************ //