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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::PointHit
29
30Description
31 Describes the interaction of a face and a point.
32 It carries the info of a successful hit and (if successful)
33 the interaction point.
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef PointHit_H
38#define PointHit_H
39
40#include "bool.H"
41#include "point.H"
42#include "token.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48
49/*---------------------------------------------------------------------------*\
50 Class PointHit Declaration
51\*---------------------------------------------------------------------------*/
52
53template<class PointType>
54class PointHit
55{
56 // Private Data
57
58 //- Point of hit; for miss holds best estimate outside the object
59 PointType point_;
60
61 //- Distance to hit point
62 scalar distance_;
63
64 //- Hit success
65 bool hit_;
66
67 //- Eligible miss
68 bool eligibleMiss_;
69
70
71public:
72
73 // Public Typedefs
74
75 //- The point type
76 typedef PointType point_type;
77
78
79 // Constructors
80
81 //- Default construct. A zero point with a large distance,
82 //- no hit, no eligible misses
83 PointHit()
84 :
85 point_(Zero),
86 distance_(GREAT),
87 hit_(false),
88 eligibleMiss_(false)
89 {}
90
91 //- Construct from point with a large distance,
92 //- no hit, no eligible misses
93 explicit PointHit(const point_type& p)
94 :
95 point_(p),
96 distance_(GREAT),
97 hit_(false),
98 eligibleMiss_(false)
99 {}
100
101 //- Construct from components
103 (
104 const bool hit,
105 const point_type& p,
106 const scalar dist,
107 const bool eligibleMiss = false
108 )
109 :
110 point_(p),
111 distance_(dist),
112 hit_(hit),
113 eligibleMiss_(eligibleMiss)
114 {}
115
116
117 // Member Functions
118
119 // Access
120
121 //- Is there a hit
122 bool hit() const noexcept
123 {
124 return hit_;
125 }
126
127 //- Is this an eligible miss
128 bool eligibleMiss() const noexcept
129 {
130 return eligibleMiss_;
131 }
132
133 //- Return the point, no checks
134 const point_type& point() const noexcept
135 {
136 return point_;
137 }
138
139 //- Return distance to hit
140 scalar distance() const noexcept
141 {
142 return distance_;
143 }
144
145 //- Return the hit point. Fatal if not hit.
146 const point_type& hitPoint() const
147 {
148 if (!hit_)
149 {
151 << "Requested a hit point, but it was not hit"
152 << abort(FatalError);
153 }
154
155 return point_;
156 }
157
158 //- Return the miss point. Fatal if hit.
159 const point_type& missPoint() const
160 {
161 if (hit_)
162 {
164 << "Requested a miss point, but it was hit"
165 << abort(FatalError);
166 }
167
168 return point_;
169 }
170
171 //- The point, no checks
172 // \deprecated(2020-10) use point()
173 const point_type& rawPoint() const noexcept
174 {
175 return point_;
176 }
177
178
179 // Edit
180
181 //- Set the hit status \em on
182 void setHit() noexcept
183 {
184 hit_ = true;
185 eligibleMiss_ = false;
186 }
187
188 //- Set the hit status \em off and set the eligible miss status
189 void setMiss(const bool eligible) noexcept
190 {
191 hit_ = false;
192 eligibleMiss_ = eligible;
193 }
194
195 //- Set the point
196 void setPoint(const point_type& p)
197 {
198 point_ = p;
199 }
200
201 //- Set the distance
202 void setDistance(const scalar d) noexcept
203 {
204 distance_ = d;
205 }
206
207
208 // Member Operators
209
210 //- Distance comparision operator, for sorting
211 template<class AnyPointType>
212 bool operator<(const PointHit<AnyPointType>& rhs) const noexcept
213 {
214 return distance_ < rhs.distance_;
215 }
216};
217
218
219// Ostream Operator
220
221template<class PointType>
222inline Ostream& operator<<(Ostream& os, const PointHit<PointType>& pHit)
223{
224 os << pHit.hit() << token::SPACE
225 << pHit.point() << token::SPACE
226 << pHit.distance() << token::SPACE
227 << pHit.eligibleMiss();
228
229 return os;
230}
231
232
233// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234
235} // End namespace Foam
236
237// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238
239#endif
240
241// ************************************************************************* //
System bool.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:54
void setHit() noexcept
Set the hit status on.
Definition: PointHit.H:181
PointType point_type
The point type.
Definition: PointHit.H:75
const point_type & rawPoint() const noexcept
The point, no checks.
Definition: PointHit.H:172
bool eligibleMiss() const noexcept
Is this an eligible miss.
Definition: PointHit.H:127
bool operator<(const PointHit< AnyPointType > &rhs) const noexcept
Distance comparision operator, for sorting.
Definition: PointHit.H:211
scalar distance() const noexcept
Return distance to hit.
Definition: PointHit.H:139
PointHit(const bool hit, const point_type &p, const scalar dist, const bool eligibleMiss=false)
Construct from components.
Definition: PointHit.H:102
void setPoint(const point_type &p)
Set the point.
Definition: PointHit.H:195
void setDistance(const scalar d) noexcept
Set the distance.
Definition: PointHit.H:201
const point_type & missPoint() const
Return the miss point. Fatal if hit.
Definition: PointHit.H:158
bool hit() const noexcept
Is there a hit.
Definition: PointHit.H:121
PointHit(const point_type &p)
Definition: PointHit.H:92
void setMiss(const bool eligible) noexcept
Set the hit status off and set the eligible miss status.
Definition: PointHit.H:188
const point_type & point() const noexcept
Return the point, no checks.
Definition: PointHit.H:133
const point_type & hitPoint() const
Return the hit point. Fatal if not hit.
Definition: PointHit.H:145
@ SPACE
Space [isspace].
Definition: token.H:125
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const direction noexcept
Definition: Scalar.H:223
error FatalError