Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABoundingBox.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// Author: HoangTRAN, 20/2/2019
28
29#ifndef G4DNABoundingBox_hh
30#define G4DNABoundingBox_hh 1
31
32#include "globals.hh"
33#include "G4ThreeVector.hh"
35{
36 public:
37 G4DNABoundingBox() = default;
38 ~G4DNABoundingBox() = default;
39 template <typename Iterator>
40 G4DNABoundingBox(Iterator begin, Iterator end);
41 G4DNABoundingBox(const std::initializer_list<G4double>& l);
42 G4DNABoundingBox(G4DNABoundingBox&& rhs) noexcept;
46 G4bool contains(const G4DNABoundingBox& other) const;
47 G4bool contains(const G4ThreeVector& point) const;
48 G4bool overlap(const G4DNABoundingBox& other, G4DNABoundingBox* out) const;
49 G4bool overlap(const G4ThreeVector& query, const G4double& radius) const;
50 G4bool contains(const G4ThreeVector& query, const G4ThreeVector& Point,
51 const G4double& Radius) const;
52 G4bool contains(const G4ThreeVector& query, const G4double& radius) const;
53 G4double Volume() const;
54 std::array<G4DNABoundingBox, 8> partition() const;
55 friend std::ostream& operator<<(std::ostream& s, const G4DNABoundingBox& rhs);
56 G4bool operator==(const G4DNABoundingBox& rhs) const;
57 G4bool operator!=(const G4DNABoundingBox& rhs) const;
58 void resize(G4ThreeVector pics[8]);
59 G4DNABoundingBox translate(const G4ThreeVector& trans) const;
60 inline G4double halfSideLengthInX() const;
61 inline G4double halfSideLengthInY() const;
62 inline G4double halfSideLengthInZ() const;
63 inline G4double Getxhi() const;
64 inline G4double Getyhi() const;
65 inline G4double Getzhi() const;
66 inline G4double Getxlo() const;
67 inline G4double Getylo() const;
68 inline G4double Getzlo() const;
69 inline G4ThreeVector middlePoint() const;
70
71 private:
72 G4double fxhi, fxlo, fyhi, fylo, fzhi, fzlo;
73};
74
76 -std::numeric_limits<G4double>::max(), std::numeric_limits<G4double>::max(),
77 -std::numeric_limits<G4double>::max(), std::numeric_limits<G4double>::max(),
78 -std::numeric_limits<G4double>::max(), std::numeric_limits<G4double>::max()
79};
81 G4DNABoundingBox{ std::numeric_limits<G4double>::quiet_NaN(),
82 std::numeric_limits<G4double>::quiet_NaN(),
83 std::numeric_limits<G4double>::quiet_NaN(),
84 std::numeric_limits<G4double>::quiet_NaN(),
85 std::numeric_limits<G4double>::quiet_NaN(),
86 std::numeric_limits<G4double>::quiet_NaN() };
87
89{
90 G4ThreeVector mid((fxhi + fxlo) / 2.0, (fyhi + fylo) / 2.0,
91 (fzhi + fzlo) / 2.0);
92 return mid;
93}
94
96{
97 return std::abs(fxhi - fxlo) * 0.5;
98}
99
101{
102 return std::abs(fyhi - fylo) * 0.5;
103}
104
106{
107 return std::abs(fzhi - fzlo) * 0.5;
108}
109
110template <typename Iterator>
111G4DNABoundingBox::G4DNABoundingBox(Iterator begin, Iterator end)
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}
146
147inline G4double G4DNABoundingBox::Getxhi() const { return fxhi; }
148inline G4double G4DNABoundingBox::Getyhi() const { return fyhi; }
149inline G4double G4DNABoundingBox::Getzhi() const { return fzhi; }
150inline G4double G4DNABoundingBox::Getxlo() const { return fxlo; }
151inline G4double G4DNABoundingBox::Getylo() const { return fylo; }
152inline G4double G4DNABoundingBox::Getzlo() const { return fzlo; }
153#endif
const G4DNABoundingBox initial
const G4DNABoundingBox invalid
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
double z() const
double x() const
double y() const
~G4DNABoundingBox()=default
void resize(G4ThreeVector pics[8])
G4bool operator!=(const G4DNABoundingBox &rhs) const
G4double Getxlo() const
G4DNABoundingBox(const std::initializer_list< G4double > &l)
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
G4DNABoundingBox(const G4DNABoundingBox &)=default
std::array< G4DNABoundingBox, 8 > partition() const
friend std::ostream & operator<<(std::ostream &s, const G4DNABoundingBox &rhs)
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