File indexing completed on 2025-01-18 09:58:07
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #ifndef G4DNABoundingBox_hh
0030 #define G4DNABoundingBox_hh 1
0031
0032 #include "globals.hh"
0033 #include "G4ThreeVector.hh"
0034 class G4DNABoundingBox
0035 {
0036 public:
0037 G4DNABoundingBox() = default;
0038 ~G4DNABoundingBox() = default;
0039 template <typename Iterator>
0040 G4DNABoundingBox(Iterator begin, Iterator end);
0041 G4DNABoundingBox(const std::initializer_list<G4double>& l);
0042 G4DNABoundingBox(G4DNABoundingBox&& rhs) noexcept;
0043 G4DNABoundingBox(const G4DNABoundingBox&) = default;
0044 G4DNABoundingBox& operator=(const G4DNABoundingBox&) = default;
0045 G4DNABoundingBox& operator=(G4DNABoundingBox&& hrs) noexcept;
0046 G4bool contains(const G4DNABoundingBox& other) const;
0047 G4bool contains(const G4ThreeVector& point) const;
0048 G4bool overlap(const G4DNABoundingBox& other, G4DNABoundingBox* out) const;
0049 G4bool overlap(const G4ThreeVector& query, const G4double& radius) const;
0050 G4bool contains(const G4ThreeVector& query, const G4ThreeVector& Point,
0051 const G4double& Radius) const;
0052 G4bool contains(const G4ThreeVector& query, const G4double& radius) const;
0053 G4double Volume() const;
0054 std::array<G4DNABoundingBox, 8> partition() const;
0055 friend std::ostream& operator<<(std::ostream& s, const G4DNABoundingBox& rhs);
0056 G4bool operator==(const G4DNABoundingBox& rhs) const;
0057 G4bool operator!=(const G4DNABoundingBox& rhs) const;
0058 void resize(G4ThreeVector pics[8]);
0059 G4DNABoundingBox translate(const G4ThreeVector& trans) const;
0060 inline G4double halfSideLengthInX() const;
0061 inline G4double halfSideLengthInY() const;
0062 inline G4double halfSideLengthInZ() const;
0063 inline G4double Getxhi() const;
0064 inline G4double Getyhi() const;
0065 inline G4double Getzhi() const;
0066 inline G4double Getxlo() const;
0067 inline G4double Getylo() const;
0068 inline G4double Getzlo() const;
0069 inline G4ThreeVector middlePoint() const;
0070
0071 private:
0072 G4double fxhi, fxlo, fyhi, fylo, fzhi, fzlo;
0073 };
0074
0075 const G4DNABoundingBox initial = G4DNABoundingBox{
0076 -std::numeric_limits<G4double>::max(), std::numeric_limits<G4double>::max(),
0077 -std::numeric_limits<G4double>::max(), std::numeric_limits<G4double>::max(),
0078 -std::numeric_limits<G4double>::max(), std::numeric_limits<G4double>::max()
0079 };
0080 const G4DNABoundingBox invalid =
0081 G4DNABoundingBox{ std::numeric_limits<G4double>::quiet_NaN(),
0082 std::numeric_limits<G4double>::quiet_NaN(),
0083 std::numeric_limits<G4double>::quiet_NaN(),
0084 std::numeric_limits<G4double>::quiet_NaN(),
0085 std::numeric_limits<G4double>::quiet_NaN(),
0086 std::numeric_limits<G4double>::quiet_NaN() };
0087
0088 inline G4ThreeVector G4DNABoundingBox::middlePoint() const
0089 {
0090 G4ThreeVector mid((fxhi + fxlo) / 2.0, (fyhi + fylo) / 2.0,
0091 (fzhi + fzlo) / 2.0);
0092 return mid;
0093 }
0094
0095 inline G4double G4DNABoundingBox::halfSideLengthInX() const
0096 {
0097 return std::abs(fxhi - fxlo) * 0.5;
0098 }
0099
0100 inline G4double G4DNABoundingBox::halfSideLengthInY() const
0101 {
0102 return std::abs(fyhi - fylo) * 0.5;
0103 }
0104
0105 inline G4double G4DNABoundingBox::halfSideLengthInZ() const
0106 {
0107 return std::abs(fzhi - fzlo) * 0.5;
0108 }
0109
0110 template <typename Iterator>
0111 G4DNABoundingBox::G4DNABoundingBox(Iterator begin, Iterator end)
0112 : G4DNABoundingBox(initial)
0113 {
0114 for(; begin != end; ++begin)
0115 {
0116 const G4ThreeVector& point = (*begin);
0117
0118 if(point.x() < fxlo)
0119 {
0120 fxlo = point.x();
0121 }
0122 if(point.x() > fxhi)
0123 {
0124 fxhi = point.x();
0125 }
0126
0127 if(point.y() < fylo)
0128 {
0129 fylo = point.y();
0130 }
0131 if(point.y() > fyhi)
0132 {
0133 fyhi = point.y();
0134 }
0135
0136 if(point.z() < fzlo)
0137 {
0138 fzlo = point.z();
0139 }
0140 if(point.z() > fzhi)
0141 {
0142 fzhi = point.z();
0143 }
0144 }
0145 }
0146
0147 inline G4double G4DNABoundingBox::Getxhi() const { return fxhi; }
0148 inline G4double G4DNABoundingBox::Getyhi() const { return fyhi; }
0149 inline G4double G4DNABoundingBox::Getzhi() const { return fzhi; }
0150 inline G4double G4DNABoundingBox::Getxlo() const { return fxlo; }
0151 inline G4double G4DNABoundingBox::Getylo() const { return fylo; }
0152 inline G4double G4DNABoundingBox::Getzlo() const { return fzlo; }
0153 #endif