Geant4-11
Public Member Functions | Private Attributes | Friends
G4DNABoundingBox Class Reference

#include <G4DNABoundingBox.hh>

Public Member Functions

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

Private Attributes

G4double fxhi
 
G4double fxlo
 
G4double fyhi
 
G4double fylo
 
G4double fzhi
 
G4double fzlo
 

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

References fxhi, fxlo, fyhi, fylo, fzhi, fzlo, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ 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 127 of file G4DNABoundingBox.cc.

128{
129 return fxlo <= other.fxlo && fxhi >= other.fxhi && fylo <= other.fylo &&
130 fyhi >= other.fyhi && fzlo <= other.fzlo && fzhi >= other.fzhi;
131}

References fxhi, fxlo, fyhi, fylo, fzhi, and fzlo.

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

◆ contains() [2/4]

G4bool G4DNABoundingBox::contains ( const G4ThreeVector point) const

Definition at line 135 of file G4DNABoundingBox.cc.

136{
137 return fxlo <= point.x() && fxhi >= point.x() && fylo <= point.y() &&
138 fyhi >= point.y() && fzlo <= point.z() && fzhi >= point.z();
139}

References fxhi, fxlo, fyhi, fylo, fzhi, fzlo, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ contains() [3/4]

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

Definition at line 215 of file G4DNABoundingBox.cc.

217{
218 G4ThreeVector temp = query - this->middlePoint();
219 G4double x = std::abs(temp.x()) + this->halfSideLengthInX();
220 G4double y = std::abs(temp.y()) + this->halfSideLengthInY();
221 G4double z = std::abs(temp.z()) + this->halfSideLengthInZ();
222 G4double norm = std::sqrt(x * x + y * y + z * z);
223 return (norm < Radius);
224}
double G4double
Definition: G4Types.hh:83
G4ThreeVector middlePoint() const
G4double halfSideLengthInY() const
G4double halfSideLengthInZ() const
G4double halfSideLengthInX() const

References halfSideLengthInX(), halfSideLengthInY(), halfSideLengthInZ(), middlePoint(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ contains() [4/4]

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

Definition at line 227 of file G4DNABoundingBox.cc.

230{
231 return (((query - Point).mag()) < Radius);
232}

◆ Getxhi()

G4double G4DNABoundingBox::Getxhi ( ) const
inline

◆ Getxlo()

G4double G4DNABoundingBox::Getxlo ( ) const
inline

◆ Getyhi()

G4double G4DNABoundingBox::Getyhi ( ) const
inline

◆ Getylo()

G4double G4DNABoundingBox::Getylo ( ) const
inline

◆ Getzhi()

G4double G4DNABoundingBox::Getzhi ( ) const
inline

◆ 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}

References fxhi, and fxlo.

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}

References fyhi, and fylo.

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}

References fzhi, and fzlo.

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}

References fxhi, fxlo, fyhi, fylo, fzhi, and fzlo.

Referenced by contains(), and overlap().

◆ operator!=()

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

Definition at line 269 of file G4DNABoundingBox.cc.

270{
271 return !operator==(rhs);
272}
G4bool operator==(const G4DNABoundingBox &rhs) const

References operator==().

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

Definition at line 52 of file G4DNABoundingBox.cc.

53{
54 fxhi = rhs.fxhi;
55 fxlo = rhs.fxlo;
56 fyhi = rhs.fyhi;
57 fylo = rhs.fylo;
58 fzhi = rhs.fzhi;
59 fzlo = rhs.fzlo;
60 return *this;
61}

References fxhi.

◆ operator==()

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

Definition at line 259 of file G4DNABoundingBox.cc.

260{
261 return (fxhi == rhs.fxhi && fxlo == rhs.fxlo && fyhi == rhs.fyhi &&
262 fylo == rhs.fylo && fzhi == rhs.fzhi && fzlo == rhs.fzlo) ||
263 (std::isnan(fxhi) && std::isnan(rhs.fxhi) && std::isnan(fxlo) &&
264 std::isnan(rhs.fxlo) && std::isnan(fyhi) && std::isnan(rhs.fyhi) &&
265 std::isnan(fylo) && std::isnan(rhs.fylo) && std::isnan(fzhi) &&
266 std::isnan(rhs.fzhi) && std::isnan(fzlo) && std::isnan(rhs.fzlo));
267}

References fxhi, fxlo, fyhi, fylo, fzhi, and fzlo.

Referenced by operator!=().

◆ overlap() [1/2]

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

Definition at line 143 of file G4DNABoundingBox.cc.

145{
146 if(contains(other))
147 {
148 *out = other;
149 return true;
150 }
151 else if(other.contains(*this))
152 {
153 *out = *this;
154 return true;
155 }
156
157 // Check if there is no intersection
158 if(fxhi < other.fxlo || fxlo > other.fxhi || fyhi < other.fylo ||
159 fylo > other.fyhi || fzhi < other.fzlo || fzlo > other.fzhi)
160 {
161 *out = invalid;
162 return false;
163 }
164
165 G4double upperX = std::min(fxhi, other.fxhi);
166 G4double upperY = std::min(fyhi, other.fyhi);
167 G4double upperZ = std::min(fzhi, other.fzhi);
168
169 G4double lowerX = std::max(fxlo, other.fxlo);
170 G4double lowerY = std::max(fylo, other.fylo);
171 G4double lowerZ = std::max(fzlo, other.fzlo);
172
173 *out = G4DNABoundingBox{ upperX, lowerX, upperY, lowerY, upperZ, lowerZ };
174 return true;
175}
const G4DNABoundingBox invalid
G4bool contains(const G4DNABoundingBox &other) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References contains(), fxhi, fxlo, fyhi, fylo, fzhi, fzlo, invalid, G4INCL::Math::max(), and G4INCL::Math::min().

◆ overlap() [2/2]

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

Definition at line 179 of file G4DNABoundingBox.cc.

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

References halfSideLengthInX(), halfSideLengthInY(), halfSideLengthInZ(), G4INCL::Math::max(), middlePoint(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ partition()

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

Definition at line 236 of file G4DNABoundingBox.cc.

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

References fxhi, fxlo, fyhi, fylo, fzhi, and fzlo.

◆ resize()

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

Definition at line 74 of file G4DNABoundingBox.cc.

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

References fxhi, fxlo, fyhi, fylo, fzhi, fzlo, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ translate()

G4DNABoundingBox G4DNABoundingBox::translate ( const G4ThreeVector trans) const

Definition at line 110 of file G4DNABoundingBox.cc.

111{
112 G4DNABoundingBox output;
113
114 output.fxhi = fxhi + trans.x();
115 output.fxlo = fxlo + trans.x();
116
117 output.fyhi = fyhi + trans.y();
118 output.fylo = fylo + trans.x();
119
120 output.fzhi = fzhi + trans.z();
121 output.fzlo = fzlo + trans.z();
122
123 return output;
124}

References fxhi, fxlo, fyhi, fylo, fzhi, fzlo, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ Volume()

G4double G4DNABoundingBox::Volume ( ) const

Definition at line 63 of file G4DNABoundingBox.cc.

64{
65 auto LengthY = this->Getyhi() - this->Getylo();
66 auto LengthX = this->Getxhi() - this->Getxlo();
67 auto LengthZ = this->Getzhi() - this->Getzlo();
68 G4double V = LengthY * LengthX * LengthZ;
69 return V;
70}
G4double Getxlo() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4double Getzlo() const
G4double Getzhi() const

References Getxhi(), Getxlo(), Getyhi(), Getylo(), Getzhi(), and Getzlo().

Referenced by G4DNAScavengerMaterial::Dump(), G4DNAGillespieDirectMethod::FindScavenging(), and G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength().

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 275 of file G4DNABoundingBox.cc.

276{
277 stream << "{" << G4BestUnit(rhs.fxhi, "Length") << ", "
278 << G4BestUnit(rhs.fxlo, "Length") << ", "
279 << G4BestUnit(rhs.fyhi, "Length") << ", "
280 << G4BestUnit(rhs.fylo, "Length") << ", "
281 << G4BestUnit(rhs.fzhi, "Length") << ", "
282 << G4BestUnit(rhs.fzlo, "Length") << ", "
283 << "}";
284 return stream;
285}
#define G4BestUnit(a, b)

Field Documentation

◆ fxhi

G4double G4DNABoundingBox::fxhi
private

◆ fxlo

G4double G4DNABoundingBox::fxlo
private

◆ fyhi

G4double G4DNABoundingBox::fyhi
private

◆ fylo

G4double G4DNABoundingBox::fylo
private

◆ fzhi

G4double G4DNABoundingBox::fzhi
private

◆ fzlo

G4double G4DNABoundingBox::fzlo
private

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