Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABoundingBox Class Reference

#include <G4DNABoundingBox.hh>

Public Member Functions

 G4DNABoundingBox ()=default
 
 ~G4DNABoundingBox ()=default
 
template<typename Iterator >
 G4DNABoundingBox (Iterator begin, Iterator end)
 
 G4DNABoundingBox (const std::initializer_list< G4double > &l)
 
 G4DNABoundingBox (G4DNABoundingBox &&rhs) noexcept
 
 G4DNABoundingBox (const G4DNABoundingBox &)=default
 
G4DNABoundingBoxoperator= (const G4DNABoundingBox &)=default
 
G4DNABoundingBoxoperator= (G4DNABoundingBox &&hrs) noexcept
 
G4bool contains (const G4DNABoundingBox &other) const
 
G4bool contains (const G4ThreeVector &point) const
 
G4bool overlap (const G4DNABoundingBox &other, G4DNABoundingBox *out) const
 
G4bool overlap (const G4ThreeVector &query, const G4double &radius) const
 
G4bool contains (const G4ThreeVector &query, const G4ThreeVector &Point, const G4double &Radius) const
 
G4bool contains (const G4ThreeVector &query, const G4double &radius) const
 
G4double Volume () const
 
std::array< G4DNABoundingBox, 8 > partition () const
 
G4bool operator== (const G4DNABoundingBox &rhs) const
 
G4bool operator!= (const G4DNABoundingBox &rhs) const
 
void resize (G4ThreeVector pics[8])
 
G4DNABoundingBox translate (const G4ThreeVector &trans) const
 
G4double halfSideLengthInX () const
 
G4double halfSideLengthInY () const
 
G4double halfSideLengthInZ () const
 
G4double Getxhi () const
 
G4double Getyhi () const
 
G4double Getzhi () const
 
G4double Getxlo () const
 
G4double Getylo () const
 
G4double Getzlo () const
 
G4ThreeVector middlePoint () const
 

Friends

std::ostream & operator<< (std::ostream &s, const G4DNABoundingBox &rhs)
 

Detailed Description

Definition at line 34 of file G4DNABoundingBox.hh.

Constructor & Destructor Documentation

◆ G4DNABoundingBox() [1/5]

G4DNABoundingBox::G4DNABoundingBox ( )
default

◆ ~G4DNABoundingBox()

G4DNABoundingBox::~G4DNABoundingBox ( )
default

◆ G4DNABoundingBox() [2/5]

template<typename Iterator >
G4DNABoundingBox::G4DNABoundingBox ( Iterator begin,
Iterator end )

Definition at line 111 of file G4DNABoundingBox.hh.

113{
114 for(; begin != end; ++begin)
115 {
116 const G4ThreeVector& point = (*begin);
117
118 if(point.x() < fxlo)
119 {
120 fxlo = point.x();
121 }
122 if(point.x() > fxhi)
123 {
124 fxhi = point.x();
125 }
126
127 if(point.y() < fylo)
128 {
129 fylo = point.y();
130 }
131 if(point.y() > fyhi)
132 {
133 fyhi = point.y();
134 }
135
136 if(point.z() < fzlo)
137 {
138 fzlo = point.z();
139 }
140 if(point.z() > fzhi)
141 {
142 fzhi = point.z();
143 }
144 }
145}
const G4DNABoundingBox initial
double z() const
double x() const
double y() const
G4DNABoundingBox()=default

◆ G4DNABoundingBox() [3/5]

G4DNABoundingBox::G4DNABoundingBox ( const std::initializer_list< G4double > & l)

◆ G4DNABoundingBox() [4/5]

G4DNABoundingBox::G4DNABoundingBox ( G4DNABoundingBox && rhs)
noexcept

Definition at line 34 of file G4DNABoundingBox.cc.

35{
36 fxhi = rhs.fxhi;
37 fxlo = rhs.fxlo;
38 fyhi = rhs.fyhi;
39 fylo = rhs.fylo;
40 fzhi = rhs.fzhi;
41 fzlo = rhs.fzlo;
42}

◆ G4DNABoundingBox() [5/5]

G4DNABoundingBox::G4DNABoundingBox ( const G4DNABoundingBox & )
default

Member Function Documentation

◆ contains() [1/4]

G4bool G4DNABoundingBox::contains ( const G4DNABoundingBox & other) const

Definition at line 133 of file G4DNABoundingBox.cc.

134{
135 return fxlo <= other.fxlo && fxhi >= other.fxhi && fylo <= other.fylo &&
136 fyhi >= other.fyhi && fzlo <= other.fzlo && fzhi >= other.fzhi;
137}

Referenced by G4DNAMesh::GetIndex(), and overlap().

◆ contains() [2/4]

G4bool G4DNABoundingBox::contains ( const G4ThreeVector & point) const

Definition at line 141 of file G4DNABoundingBox.cc.

142{
143 return fxlo <= point.x() && fxhi >= point.x() && fylo <= point.y() &&
144 fyhi >= point.y() && fzlo <= point.z() && fzhi >= point.z();
145}

◆ contains() [3/4]

G4bool G4DNABoundingBox::contains ( const G4ThreeVector & query,
const G4double & radius ) const

Definition at line 221 of file G4DNABoundingBox.cc.

223{
224 G4ThreeVector temp = query - this->middlePoint();
225 G4double x = std::abs(temp.x()) + this->halfSideLengthInX();
226 G4double y = std::abs(temp.y()) + this->halfSideLengthInY();
227 G4double z = std::abs(temp.z()) + this->halfSideLengthInZ();
228 G4double norm = std::sqrt(x * x + y * y + z * z);
229 return (norm < Radius);
230}
double G4double
Definition G4Types.hh:83
G4ThreeVector middlePoint() const
G4double halfSideLengthInY() const
G4double halfSideLengthInZ() const
G4double halfSideLengthInX() const

◆ contains() [4/4]

G4bool G4DNABoundingBox::contains ( const G4ThreeVector & query,
const G4ThreeVector & Point,
const G4double & Radius ) const

Definition at line 233 of file G4DNABoundingBox.cc.

236{
237 return (((query - Point).mag()) < Radius);
238}

◆ Getxhi()

G4double G4DNABoundingBox::Getxhi ( ) const
inline

Definition at line 147 of file G4DNABoundingBox.hh.

147{ return fxhi; }

Referenced by G4DNAMesh::ConvertIndex(), G4DNAMesh::FindNeighboringVoxels(), and Volume().

◆ Getxlo()

G4double G4DNABoundingBox::Getxlo ( ) const
inline

◆ Getyhi()

G4double G4DNABoundingBox::Getyhi ( ) const
inline

Definition at line 148 of file G4DNABoundingBox.hh.

148{ return fyhi; }

Referenced by G4DNAMesh::ConvertIndex(), G4DNAMesh::FindNeighboringVoxels(), and Volume().

◆ Getylo()

G4double G4DNABoundingBox::Getylo ( ) const
inline

◆ Getzhi()

G4double G4DNABoundingBox::Getzhi ( ) const
inline

Definition at line 149 of file G4DNABoundingBox.hh.

149{ return fzhi; }

Referenced by G4DNAMesh::ConvertIndex(), G4DNAMesh::FindNeighboringVoxels(), and Volume().

◆ Getzlo()

G4double G4DNABoundingBox::Getzlo ( ) const
inline

◆ halfSideLengthInX()

G4double G4DNABoundingBox::halfSideLengthInX ( ) const
inline

Definition at line 95 of file G4DNABoundingBox.hh.

96{
97 return std::abs(fxhi - fxlo) * 0.5;
98}

Referenced by contains(), and overlap().

◆ halfSideLengthInY()

G4double G4DNABoundingBox::halfSideLengthInY ( ) const
inline

Definition at line 100 of file G4DNABoundingBox.hh.

101{
102 return std::abs(fyhi - fylo) * 0.5;
103}

Referenced by contains(), and overlap().

◆ halfSideLengthInZ()

G4double G4DNABoundingBox::halfSideLengthInZ ( ) const
inline

Definition at line 105 of file G4DNABoundingBox.hh.

106{
107 return std::abs(fzhi - fzlo) * 0.5;
108}

Referenced by contains(), and overlap().

◆ middlePoint()

G4ThreeVector G4DNABoundingBox::middlePoint ( ) const
inline

Definition at line 88 of file G4DNABoundingBox.hh.

89{
90 G4ThreeVector mid((fxhi + fxlo) / 2.0, (fyhi + fylo) / 2.0,
91 (fzhi + fzlo) / 2.0);
92 return mid;
93}

Referenced by contains(), and overlap().

◆ operator!=()

G4bool G4DNABoundingBox::operator!= ( const G4DNABoundingBox & rhs) const

Definition at line 275 of file G4DNABoundingBox.cc.

276{
277 return !operator==(rhs);
278}
G4bool operator==(const G4DNABoundingBox &rhs) const

◆ operator=() [1/2]

G4DNABoundingBox & G4DNABoundingBox::operator= ( const G4DNABoundingBox & )
default

◆ operator=() [2/2]

G4DNABoundingBox & G4DNABoundingBox::operator= ( G4DNABoundingBox && hrs)
noexcept

Definition at line 58 of file G4DNABoundingBox.cc.

59{
60 fxhi = rhs.fxhi;
61 fxlo = rhs.fxlo;
62 fyhi = rhs.fyhi;
63 fylo = rhs.fylo;
64 fzhi = rhs.fzhi;
65 fzlo = rhs.fzlo;
66 return *this;
67}

◆ operator==()

G4bool G4DNABoundingBox::operator== ( const G4DNABoundingBox & rhs) const

Definition at line 265 of file G4DNABoundingBox.cc.

266{
267 return (fxhi == rhs.fxhi && fxlo == rhs.fxlo && fyhi == rhs.fyhi &&
268 fylo == rhs.fylo && fzhi == rhs.fzhi && fzlo == rhs.fzlo) ||
269 (std::isnan(fxhi) && std::isnan(rhs.fxhi) && std::isnan(fxlo) &&
270 std::isnan(rhs.fxlo) && std::isnan(fyhi) && std::isnan(rhs.fyhi) &&
271 std::isnan(fylo) && std::isnan(rhs.fylo) && std::isnan(fzhi) &&
272 std::isnan(rhs.fzhi) && std::isnan(fzlo) && std::isnan(rhs.fzlo));
273}

Referenced by operator!=().

◆ overlap() [1/2]

G4bool G4DNABoundingBox::overlap ( const G4DNABoundingBox & other,
G4DNABoundingBox * out ) const

Definition at line 149 of file G4DNABoundingBox.cc.

151{
152 if(contains(other))
153 {
154 *out = other;
155 return true;
156 }
157 if(other.contains(*this))
158 {
159 *out = *this;
160 return true;
161 }
162
163 // Check if there is no intersection
164 if(fxhi < other.fxlo || fxlo > other.fxhi || fyhi < other.fylo ||
165 fylo > other.fyhi || fzhi < other.fzlo || fzlo > other.fzhi)
166 {
167 *out = invalid;
168 return false;
169 }
170
171 G4double upperX = std::min(fxhi, other.fxhi);
172 G4double upperY = std::min(fyhi, other.fyhi);
173 G4double upperZ = std::min(fzhi, other.fzhi);
174
175 G4double lowerX = std::max(fxlo, other.fxlo);
176 G4double lowerY = std::max(fylo, other.fylo);
177 G4double lowerZ = std::max(fzlo, other.fzlo);
178
179 *out = G4DNABoundingBox{ upperX, lowerX, upperY, lowerY, upperZ, lowerZ };
180 return true;
181}
const G4DNABoundingBox invalid
G4bool contains(const G4DNABoundingBox &other) const

◆ overlap() [2/2]

G4bool G4DNABoundingBox::overlap ( const G4ThreeVector & query,
const G4double & radius ) const

Definition at line 185 of file G4DNABoundingBox.cc.

187{
188 G4double x = query.x() - this->middlePoint().x();
189 G4double y = query.y() - this->middlePoint().y();
190 G4double z = query.z() - this->middlePoint().z();
191
192 x = std::abs(x);
193 y = std::abs(y);
194 z = std::abs(z);
195
196 if((x > (Radius + this->halfSideLengthInX())) ||
197 (y > (Radius + this->halfSideLengthInY())) ||
198 (z > (Radius + this->halfSideLengthInZ())))
199 {
200 return false; //
201 }
202 G4int num_less_extent = static_cast<int>(x < this->halfSideLengthInX()) +
203 static_cast<int>(y < this->halfSideLengthInY()) +
204 static_cast<int>(z < this->halfSideLengthInZ());
205
206 if(num_less_extent > 1)
207 {
208 return true;
209 }
210
211 x = std::max(x - this->halfSideLengthInX(), 0.0);
212 y = std::max(y - this->halfSideLengthInY(), 0.0);
213 z = std::max(z - this->halfSideLengthInZ(), 0.0);
214
215 G4double norm = std::sqrt(x * x + y * y + z * z);
216
217 return (norm < Radius);
218}
int G4int
Definition G4Types.hh:85

◆ partition()

array< G4DNABoundingBox, 8 > G4DNABoundingBox::partition ( ) const

Definition at line 242 of file G4DNABoundingBox.cc.

243{
244 G4double xmid = (fxhi + fxlo) / 2.;
245 G4double ymid = (fyhi + fylo) / 2.;
246 G4double zmid = (fzhi + fzlo) / 2.;
247
248 std::array<G4DNABoundingBox, 8> ret{ {
249 G4DNABoundingBox{ xmid, fxlo, ymid, fylo, zmid,
250 fzlo }, // bottom left front
251 G4DNABoundingBox{ fxhi, xmid, ymid, fylo, zmid,
252 fzlo }, // bottom right front
253 G4DNABoundingBox{ xmid, fxlo, fyhi, ymid, zmid, fzlo }, // bottom left back
254 G4DNABoundingBox{ fxhi, xmid, fyhi, ymid, zmid,
255 fzlo }, // bottom right back
256 G4DNABoundingBox{ xmid, fxlo, ymid, fylo, fzhi, zmid }, // top left front
257 G4DNABoundingBox{ fxhi, xmid, ymid, fylo, fzhi, zmid }, // top right front
258 G4DNABoundingBox{ xmid, fxlo, fyhi, ymid, fzhi, zmid }, // top left back
259 G4DNABoundingBox{ fxhi, xmid, fyhi, ymid, fzhi, zmid } // top right back
260 } };
261 return ret;
262}

◆ resize()

void G4DNABoundingBox::resize ( G4ThreeVector pics[8])

Definition at line 80 of file G4DNABoundingBox.cc.

81{
82 for(size_t i = 0; i < 8; i++)
83 {
84 const G4ThreeVector& point = pics[i];
85 if(point.x() < fxlo)
86 {
87 fxlo = point.x();
88 }
89 if(point.x() > fxhi)
90 {
91 fxhi = point.x();
92 }
93
94 if(point.y() < fylo)
95 {
96 fylo = point.y();
97 }
98 if(point.y() > fyhi)
99 {
100 fyhi = point.y();
101 }
102
103 if(point.z() < fzlo)
104 {
105 fzlo = point.z();
106 }
107 if(point.z() > fzhi)
108 {
109 fzhi = point.z();
110 }
111 }
112}

◆ translate()

G4DNABoundingBox G4DNABoundingBox::translate ( const G4ThreeVector & trans) const

Definition at line 116 of file G4DNABoundingBox.cc.

117{
118 G4DNABoundingBox output;
119
120 output.fxhi = fxhi + trans.x();
121 output.fxlo = fxlo + trans.x();
122
123 output.fyhi = fyhi + trans.y();
124 output.fylo = fylo + trans.x();
125
126 output.fzhi = fzhi + trans.z();
127 output.fzlo = fzlo + trans.z();
128
129 return output;
130}

◆ Volume()

G4double G4DNABoundingBox::Volume ( ) const

Definition at line 69 of file G4DNABoundingBox.cc.

70{
71 auto LengthY = this->Getyhi() - this->Getylo();
72 auto LengthX = this->Getxhi() - this->Getxlo();
73 auto LengthZ = this->Getzhi() - this->Getzlo();
74 G4double V = LengthY * LengthX * LengthZ;
75 return V;
76}
G4double Getxlo() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4double Getzlo() const
G4double Getzhi() const

Referenced by G4DNAScavengerMaterial::Dump(), G4DNAScavengerMaterial::GetpH(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), and G4DNAScavengerMaterial::SetpH().

Friends And Related Symbol Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream & s,
const G4DNABoundingBox & rhs )
friend

Definition at line 281 of file G4DNABoundingBox.cc.

282{
283 stream << "{" << G4BestUnit(rhs.fxhi, "Length") << ", "
284 << G4BestUnit(rhs.fxlo, "Length") << ", "
285 << G4BestUnit(rhs.fyhi, "Length") << ", "
286 << G4BestUnit(rhs.fylo, "Length") << ", "
287 << G4BestUnit(rhs.fzhi, "Length") << ", "
288 << G4BestUnit(rhs.fzlo, "Length") << ", "
289 << "}";
290 return stream;
291}
#define G4BestUnit(a, b)

The documentation for this class was generated from the following files: