Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/management/src/G4DNABoundingBox.cc

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 // * 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"
 30 using std::array;
 31 using std::initializer_list;
 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 33 
 34 G4DNABoundingBox::G4DNABoundingBox(G4DNABoundingBox&& rhs) noexcept
 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 }
 43 ////....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 44 
 45 G4DNABoundingBox::G4DNABoundingBox(const initializer_list<G4double>& l)
 46 {
 47   fxhi = 0.;
 48   fxlo = 0.;
 49   fyhi = 0.;
 50   fylo = 0.;
 51   fzhi = 0.;
 52   fzlo = 0.;
 53   std::copy(l.begin(), l.end(), &fxhi);
 54 }
 55 
 56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 57 
 58 G4DNABoundingBox& G4DNABoundingBox::operator=(G4DNABoundingBox&& rhs) noexcept
 59 {
 60   fxhi = rhs.fxhi;
 61   fxlo = rhs.fxlo;
 62   fyhi = rhs.fyhi;
 63   fylo = rhs.fylo;
 64   fzhi = rhs.fzhi;
 65   fzlo = rhs.fzlo;
 66   return *this;
 67 }
 68 
 69 G4double G4DNABoundingBox::Volume() const
 70 {
 71   auto LengthY = this->Getyhi() - this->Getylo();
 72   auto LengthX = this->Getxhi() - this->Getxlo();
 73   auto LengthZ = this->Getzhi() - this->Getzlo();
 74   G4double V   = LengthY * LengthX * LengthZ;
 75   return V;
 76 }
 77 
 78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 79 
 80 void G4DNABoundingBox::resize(G4ThreeVector pics[8])
 81 {
 82   for(size_t i = 0; i < 8; i++)
 83   {
 84     const G4ThreeVector& point = pics[i];
 85     if(point.x() < fxlo)
 86     {
 87       fxlo = point.x();
 88     }
 89     if(point.x() > fxhi)
 90     {
 91       fxhi = point.x();
 92     }
 93 
 94     if(point.y() < fylo)
 95     {
 96       fylo = point.y();
 97     }
 98     if(point.y() > fyhi)
 99     {
100       fyhi = point.y();
101     }
102 
103     if(point.z() < fzlo)
104     {
105       fzlo = point.z();
106     }
107     if(point.z() > fzhi)
108     {
109       fzhi = point.z();
110     }
111   }
112 }
113 
114 //------------------------------------------------------------------------------
115 
116 G4DNABoundingBox G4DNABoundingBox::translate(const G4ThreeVector& trans) const
117 {
118   G4DNABoundingBox output;
119 
120   output.fxhi = fxhi + trans.x();
121   output.fxlo = fxlo + trans.x();
122 
123   output.fyhi = fyhi + trans.y();
124   output.fylo = fylo + trans.x();
125 
126   output.fzhi = fzhi + trans.z();
127   output.fzlo = fzlo + trans.z();
128 
129   return output;
130 }
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132 
133 G4bool G4DNABoundingBox::contains(const G4DNABoundingBox& other) const
134 {
135   return fxlo <= other.fxlo && fxhi >= other.fxhi && fylo <= other.fylo &&
136          fyhi >= other.fyhi && fzlo <= other.fzlo && fzhi >= other.fzhi;
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
141 G4bool G4DNABoundingBox::contains(const G4ThreeVector& point) const
142 {
143   return fxlo <= point.x() && fxhi >= point.x() && fylo <= point.y() &&
144          fyhi >= point.y() && fzlo <= point.z() && fzhi >= point.z();
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
149 G4bool G4DNABoundingBox::overlap(const G4DNABoundingBox& other,
150                                  G4DNABoundingBox* out) const
151 {
152   if(contains(other))
153   {
154     *out = other;
155     return true;
156   }
157   if(other.contains(*this))
158   {
159     *out = *this;
160     return true;
161   }
162 
163   // Check if there is no intersection
164   if(fxhi < other.fxlo || fxlo > other.fxhi || fyhi < other.fylo ||
165      fylo > other.fyhi || fzhi < other.fzlo || fzlo > other.fzhi)
166   {
167     *out = invalid;
168     return false;
169   }
170 
171   G4double upperX = std::min(fxhi, other.fxhi);
172   G4double upperY = std::min(fyhi, other.fyhi);
173   G4double upperZ = std::min(fzhi, other.fzhi);
174 
175   G4double lowerX = std::max(fxlo, other.fxlo);
176   G4double lowerY = std::max(fylo, other.fylo);
177   G4double lowerZ = std::max(fzlo, other.fzlo);
178 
179   *out = G4DNABoundingBox{ upperX, lowerX, upperY, lowerY, upperZ, lowerZ };
180   return true;
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184 
185 G4bool G4DNABoundingBox::overlap(const G4ThreeVector& query,
186                                  const G4double& Radius) const
187 {
188   G4double x = query.x() - this->middlePoint().x();
189   G4double y = query.y() - this->middlePoint().y();
190   G4double z = query.z() - this->middlePoint().z();
191 
192   x = std::abs(x);
193   y = std::abs(y);
194   z = std::abs(z);
195 
196   if((x > (Radius + this->halfSideLengthInX())) ||
197      (y > (Radius + this->halfSideLengthInY())) ||
198      (z > (Radius + this->halfSideLengthInZ())))
199   {
200     return false;  //
201   }
202   G4int num_less_extent = static_cast<int>(x < this->halfSideLengthInX()) +
203                           static_cast<int>(y < this->halfSideLengthInY()) +
204                           static_cast<int>(z < this->halfSideLengthInZ());
205 
206   if(num_less_extent > 1)
207   {
208     return true;
209   }
210 
211   x = std::max(x - this->halfSideLengthInX(), 0.0);
212   y = std::max(y - this->halfSideLengthInY(), 0.0);
213   z = std::max(z - this->halfSideLengthInZ(), 0.0);
214 
215   G4double norm = std::sqrt(x * x + y * y + z * z);
216 
217   return (norm < Radius);
218 }
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220 
221 G4bool G4DNABoundingBox::contains(const G4ThreeVector& query,
222                                   const G4double& Radius) const
223 {
224   G4ThreeVector temp = query - this->middlePoint();
225   G4double x         = std::abs(temp.x()) + this->halfSideLengthInX();
226   G4double y         = std::abs(temp.y()) + this->halfSideLengthInY();
227   G4double z         = std::abs(temp.z()) + this->halfSideLengthInZ();
228   G4double norm      = std::sqrt(x * x + y * y + z * z);
229   return (norm < Radius);
230 }
231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
232 
233 G4bool G4DNABoundingBox::contains(const G4ThreeVector& query,
234                                   const G4ThreeVector& Point,
235                                   const G4double& Radius) const
236 {
237   return (((query - Point).mag()) < Radius);
238 }
239 
240 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
241 
242 array<G4DNABoundingBox, 8> G4DNABoundingBox::partition() const
243 {
244   G4double xmid = (fxhi + fxlo) / 2.;
245   G4double ymid = (fyhi + fylo) / 2.;
246   G4double zmid = (fzhi + fzlo) / 2.;
247 
248   std::array<G4DNABoundingBox, 8> ret{ {
249     G4DNABoundingBox{ xmid, fxlo, ymid, fylo, zmid,
250                       fzlo },  // bottom left front
251     G4DNABoundingBox{ fxhi, xmid, ymid, fylo, zmid,
252                       fzlo },  // bottom right front
253     G4DNABoundingBox{ xmid, fxlo, fyhi, ymid, zmid, fzlo },  // bottom left back
254     G4DNABoundingBox{ fxhi, xmid, fyhi, ymid, zmid,
255                       fzlo },  // bottom right back
256     G4DNABoundingBox{ xmid, fxlo, ymid, fylo, fzhi, zmid },  // top left front
257     G4DNABoundingBox{ fxhi, xmid, ymid, fylo, fzhi, zmid },  // top right front
258     G4DNABoundingBox{ xmid, fxlo, fyhi, ymid, fzhi, zmid },  // top left back
259     G4DNABoundingBox{ fxhi, xmid, fyhi, ymid, fzhi, zmid }   // top right back
260   } };
261   return ret;
262 }
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
264 
265 G4bool G4DNABoundingBox::operator==(const G4DNABoundingBox& rhs) const
266 {
267   return (fxhi == rhs.fxhi && fxlo == rhs.fxlo && fyhi == rhs.fyhi &&
268           fylo == rhs.fylo && fzhi == rhs.fzhi && fzlo == rhs.fzlo) ||
269          (std::isnan(fxhi) && std::isnan(rhs.fxhi) && std::isnan(fxlo) &&
270           std::isnan(rhs.fxlo) && std::isnan(fyhi) && std::isnan(rhs.fyhi) &&
271           std::isnan(fylo) && std::isnan(rhs.fylo) && std::isnan(fzhi) &&
272           std::isnan(rhs.fzhi) && std::isnan(fzlo) && std::isnan(rhs.fzlo));
273 }
274 
275 G4bool G4DNABoundingBox::operator!=(const G4DNABoundingBox& rhs) const
276 {
277   return !operator==(rhs);
278 }
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280 
281 std::ostream& operator<<(std::ostream& stream, const G4DNABoundingBox& rhs)
282 {
283   stream << "{" << G4BestUnit(rhs.fxhi, "Length") << ", "
284          << G4BestUnit(rhs.fxlo, "Length") << ", "
285          << G4BestUnit(rhs.fyhi, "Length") << ", "
286          << G4BestUnit(rhs.fylo, "Length") << ", "
287          << G4BestUnit(rhs.fzhi, "Length") << ", "
288          << G4BestUnit(rhs.fzlo, "Length") << ", "
289          << "}";
290   return stream;
291 }
292