Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/management/include/G4DNABoundingBox.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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"
 34 class G4DNABoundingBox
 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;
 43   G4DNABoundingBox(const G4DNABoundingBox&) = default;
 44   G4DNABoundingBox& operator=(const G4DNABoundingBox&) = default;
 45   G4DNABoundingBox& operator=(G4DNABoundingBox&& hrs) 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 
 75 const G4DNABoundingBox initial = G4DNABoundingBox{
 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 };
 80 const G4DNABoundingBox invalid =
 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 
 88 inline G4ThreeVector G4DNABoundingBox::middlePoint() const
 89 {
 90   G4ThreeVector mid((fxhi + fxlo) / 2.0, (fyhi + fylo) / 2.0,
 91                     (fzhi + fzlo) / 2.0);
 92   return mid;
 93 }
 94 
 95 inline G4double G4DNABoundingBox::halfSideLengthInX() const
 96 {
 97   return std::abs(fxhi - fxlo) * 0.5;
 98 }
 99 
100 inline G4double G4DNABoundingBox::halfSideLengthInY() const
101 {
102   return std::abs(fyhi - fylo) * 0.5;
103 }
104 
105 inline G4double G4DNABoundingBox::halfSideLengthInZ() const
106 {
107   return std::abs(fzhi - fzlo) * 0.5;
108 }
109 
110 template <typename Iterator>
111 G4DNABoundingBox::G4DNABoundingBox(Iterator begin, Iterator end)
112   : G4DNABoundingBox(initial)
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 
147 inline G4double G4DNABoundingBox::Getxhi() const { return fxhi; }
148 inline G4double G4DNABoundingBox::Getyhi() const { return fyhi; }
149 inline G4double G4DNABoundingBox::Getzhi() const { return fzhi; }
150 inline G4double G4DNABoundingBox::Getxlo() const { return fxlo; }
151 inline G4double G4DNABoundingBox::Getylo() const { return fylo; }
152 inline G4double G4DNABoundingBox::Getzlo() const { return fzlo; }
153 #endif
154