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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 Class
27  Foam::PointHit
28 
29 Description
30  This class describes the interaction of a face and a point. It
31  carries the info of a successful hit and (if successful), returns
32  the interaction point.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef PointHit_H
37 #define PointHit_H
38 
39 #include "bool.H"
40 #include "token.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // Forward declaration of classes
48 
49 class Ostream;
50 
51 
52 // Forward declaration of friend functions and operators
53 
54 template<class Point> class PointHit;
55 
56 template<class Point>
57 inline Ostream& operator<<(Ostream&, const PointHit<Point>&);
58 
59 
60 /*---------------------------------------------------------------------------*\
61  Class PointHit Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 template<class Point>
65 class PointHit
66 {
67  // Private data
68 
69  //- Hit success
70  bool hit_;
71 
72  //- Point of hit; for miss holds best estimate outside the object
73  Point hitPoint_;
74 
75  //- Distance to hit point
76  scalar distance_;
77 
78  //- Eligible miss
79  bool eligibleMiss_;
80 
81 
82 public:
83 
84  // Constructors
85 
86  //- Construct null
87  PointHit()
88  :
89  hit_(false),
90  hitPoint_(Zero),
91  distance_(GREAT),
92  eligibleMiss_(false)
93  {}
94 
95  //- Construct from components
97  (
98  const bool hit,
99  const Point& p,
100  const scalar dist,
101  const bool eligibleMiss
102  )
103  :
104  hit_(hit),
105  hitPoint_(p),
106  distance_(dist),
107  eligibleMiss_(eligibleMiss)
108  {}
109 
110  //- Construct from point. Hit and distance set later
111  PointHit(const Point& p)
112  :
113  hit_(false),
114  hitPoint_(p),
115  distance_(GREAT),
116  eligibleMiss_(false)
117  {}
118 
119 
120  // Member Functions
121 
122  // Access
123 
124  //- Is there a hit
125  bool hit() const
126  {
127  return hit_;
128  }
129 
130  //- Return hit point
131  const Point& hitPoint() const
132  {
133  if (!hit_)
134  {
136  << "requested a hit point for a miss"
137  << abort(FatalError);
138  }
139 
140  return hitPoint_;
141  }
142 
143  //- Return distance to hit
144  scalar distance() const
145  {
146  return distance_;
147  }
148 
149  //- Return miss point
150  const Point& missPoint() const
151  {
152  if (hit_)
153  {
155  << "requested a miss point for a hit"
156  << abort(FatalError);
157  }
158 
159  return hitPoint_;
160  }
161 
162  //- Return point with no checking
163  const Point& rawPoint() const
164  {
165  return hitPoint_;
166  }
167 
168  //- Is this an eligible miss
169  bool eligibleMiss() const
170  {
171  return eligibleMiss_;
172  }
173 
174  // Edit
175 
176  void setHit()
177  {
178  hit_ = true;
179  eligibleMiss_ = false;
180  }
181 
182  void setMiss(const bool eligible)
183  {
184  hit_ = false;
185  eligibleMiss_ = eligible;
186  }
187 
188  void setPoint(const Point& p)
189  {
190  hitPoint_ = p;
191  }
192 
193  void setDistance(const scalar d)
194  {
195  distance_ = d;
196  }
197 
198 
199  // Ostream operator
200 
201  friend Ostream& operator<< <Point>
202  (
203  Ostream& os,
204  const PointHit<Point>& b
205  );
206 };
207 
208 
209 template<class Point>
210 inline Ostream& operator<<(Ostream& os, const PointHit<Point>& b)
211 {
212  os << b.hit() << token::SPACE
213  << b.rawPoint() << token::SPACE
214  << b.distance() << token::SPACE
215  << b.eligibleMiss();
216 
217  return os;
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 } // End namespace Foam
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 #endif
228 
229 // ************************************************************************* //
token.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::PointHit::setMiss
void setMiss(const bool eligible)
Definition: PointHit.H:181
Foam::PointHit::hit
bool hit() const
Is there a hit.
Definition: PointHit.H:124
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:53
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::PointHit::PointHit
PointHit(const Point &p)
Construct from point. Hit and distance set later.
Definition: PointHit.H:110
Foam::PointHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:162
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:143
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::PointHit::eligibleMiss
bool eligibleMiss() const
Is this an eligible miss.
Definition: PointHit.H:168
Foam::PointHit::setPoint
void setPoint(const Point &p)
Definition: PointHit.H:187
bool.H
System bool.
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::PointHit::PointHit
PointHit()
Construct null.
Definition: PointHit.H:86
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::token::SPACE
Space [isspace].
Definition: token.H:112
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:53
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::PointHit::missPoint
const Point & missPoint() const
Return miss point.
Definition: PointHit.H:149
Foam::PointHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:130
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &)
Definition: boundaryPatch.C:102
Foam::PointHit::setHit
void setHit()
Definition: PointHit.H:175
Foam::PointHit::setDistance
void setDistance(const scalar d)
Definition: PointHit.H:192