Geant4-11
G4DNABoundingBox.cc
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25//
26// Author: HoangTRAN, 20/2/2019
27#include "G4DNABoundingBox.hh"
28#include <array>
29#include "G4UnitsTable.hh"
30using std::array;
31using std::initializer_list;
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
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}
44
45G4DNABoundingBox::G4DNABoundingBox(const initializer_list<G4double>& l)
46{
47 std::copy(l.begin(), l.end(), &fxhi);
48}
49
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51
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}
62
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}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
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}
107
108//------------------------------------------------------------------------------
109
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}
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126
128{
129 return fxlo <= other.fxlo && fxhi >= other.fxhi && fylo <= other.fylo &&
130 fyhi >= other.fyhi && fzlo <= other.fzlo && fzhi >= other.fzhi;
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
136{
137 return fxlo <= point.x() && fxhi >= point.x() && fylo <= point.y() &&
138 fyhi >= point.y() && fzlo <= point.z() && fzhi >= point.z();
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
144 G4DNABoundingBox* out) const
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}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
178
180 const G4double& Radius) const
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}
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214
216 const G4double& Radius) const
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}
225//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
226
228 const G4ThreeVector& Point,
229 const G4double& Radius) const
230{
231 return (((query - Point).mag()) < Radius);
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235
236array<G4DNABoundingBox, 8> G4DNABoundingBox::partition() const
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}
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
258
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}
268
270{
271 return !operator==(rhs);
272}
273//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
274
275std::ostream& operator<<(std::ostream& stream, const G4DNABoundingBox& rhs)
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}
const G4DNABoundingBox invalid
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
double y() const
void resize(G4ThreeVector pics[8])
G4bool operator!=(const G4DNABoundingBox &rhs) const
G4double Getxlo() const
G4double Volume() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4ThreeVector middlePoint() const
G4double halfSideLengthInY() const
G4DNABoundingBox translate(const G4ThreeVector &trans) const
G4DNABoundingBox()=default
G4bool operator==(const G4DNABoundingBox &rhs) const
std::array< G4DNABoundingBox, 8 > partition() const
G4bool overlap(const G4DNABoundingBox &other, G4DNABoundingBox *out) const
G4DNABoundingBox & operator=(const G4DNABoundingBox &)=default
G4bool contains(const G4DNABoundingBox &other) const
G4double halfSideLengthInZ() const
G4double halfSideLengthInX() const
G4double Getzlo() const
G4double Getzhi() const
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
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
void copy(G4double dst[], const G4double src[], size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98