ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ClusterSBPoints.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file ClusterSBPoints.cc
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 // This example is provided by the Geant4-DNA collaboration
27 // Any report or published results obtained using the Geant4-DNA software
28 // shall cite the following Geant4-DNA collaboration publication:
29 // Med. Phys. 37 (2010) 4692-4708
30 // The Geant4-DNA web site is available at http://geant4-dna.org
31 //
32 // Authors: Henri Payno and Yann Perrot
33 //
34 //
37 
38 #include "ClusterSBPoints.hh"
39 #include "G4SystemOfUnits.hh"
40 #include <iostream>
41 
42 using namespace std;
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45 
46 ClusterSBPoints::ClusterSBPoints(std::set<SBPoint*> pSBPoints):
47  fpRegisteredSBPoints()
48 {
50  std::set<SBPoint*>::iterator itPt;
51  for(itPt = pSBPoints.begin(); itPt != pSBPoints.end(); ++itPt)
52  {
53  AddSBPoint(*itPt);
54  }
55 }
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
60 {
61  Clear();
62 }
63 
65 {
66  std::set<SBPoint*>::iterator itPt;
67  for(itPt = fpRegisteredSBPoints.begin();
68  itPt != fpRegisteredSBPoints.end();
69  ++itPt)
70  {
71  (*itPt)->CleanCluster();
72  }
73  fpRegisteredSBPoints.clear();
74 }
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77 
79 {
80  assert(pSBPoint);
81  fpRegisteredSBPoints.insert(pSBPoint);
82  pSBPoint->SetCluster(this);
83 
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
90 {
91  G4double x = 0;
92  G4double y = 0;
93  G4double z = 0;
94 
95  std::set<SBPoint* >::iterator itSDSPt;
96  for(itSDSPt = fpRegisteredSBPoints.begin();
97  itSDSPt != fpRegisteredSBPoints.end();
98  ++itSDSPt)
99  {
100  x+= (*itSDSPt)->GetPosition().x();
101  y+= (*itSDSPt)->GetPosition().y();
102  z+= (*itSDSPt)->GetPosition().z();
103  }
104 
105  return G4ThreeVector(
106  x/fpRegisteredSBPoints.size(),
107  y/fpRegisteredSBPoints.size(),
108  z/fpRegisteredSBPoints.size());
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112 
114 {
115  G4double res = 0;
116  std::set<SBPoint* >::iterator itSDSPt;
117  for(itSDSPt = fpRegisteredSBPoints.begin();
118  itSDSPt != fpRegisteredSBPoints.end();
119  ++itSDSPt)
120  {
121  res += (*itSDSPt)->GetEdep();
122  }
123  return res;
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
129 {
130  fIsDoubleSB = false;
131  bool firstStrandTouch = false;
132  bool secondStrandTouch = false;
133 
134  std::set<SBPoint* >::iterator itSDSPt;
135  for(itSDSPt = fpRegisteredSBPoints.begin();
136  itSDSPt != fpRegisteredSBPoints.end();
137  ++itSDSPt)
138  {
139  // if the SDSPoint is localized on the first strand
140  if( ((*itSDSPt)->GetTouchedStrand() == 0 ) && !firstStrandTouch )
141  {
142  firstStrandTouch = true;
143  if(secondStrandTouch)
144  {
145  fIsDoubleSB = true;
146  return;
147  }
148  }
149  // if the SDSPoint is localized on the second strand
150  if( ((*itSDSPt)->GetTouchedStrand() == 1 ) && !secondStrandTouch )
151  {
152  secondStrandTouch = true;
153  if(firstStrandTouch)
154  {
155  fIsDoubleSB = true;
156  return;
157  }
158  }
159  }
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163 
164 bool AreOnTheSameCluster(const SBPoint* pPt1, const SBPoint* pPt2,
165  G4double pMinDist)
166 {
167  assert(pPt1);
168  assert(pPt2);
169 
170  G4double x1=pPt1->GetPosition().x()/nm;
171  G4double y1=pPt1->GetPosition().y()/nm;
172  G4double z1=pPt1->GetPosition().z()/nm;
173 
174  G4double x2=pPt2->GetPosition().x()/nm;
175  G4double y2=pPt2->GetPosition().y()/nm;
176  G4double z2=pPt2->GetPosition().z()/nm;
177 
178  // if the two points are closed enough
179  if(((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2))<=
180  (pMinDist/nm*pMinDist/nm))
181  {
182  return true;
183  }else
184  {
185  return false;
186  }
187 }
188 
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190 
191 void ClusterSBPoints::FindAllPointsPossible(std::vector<SBPoint*>* pPtsToCheck,
192  G4int pMinPts, G4double pMinDist)
193 {
194  assert((unsigned int)pMinPts > this->GetSize());
195  std::vector<SBPoint*>::iterator itPt = pPtsToCheck->begin();
196  while(itPt != pPtsToCheck->end() )
197  {
198 
199  // If 1- each SBpoint is part of only one cluster
200  // 2- the point isn't already in the cluster
201  // 3- the point is close enough of the barycenter
202  if( (!(*itPt)->HasCluster())
203  && (fpRegisteredSBPoints.find(*itPt) == fpRegisteredSBPoints.end())
204  && HasInBarycenter(*itPt, pMinDist)) // first version used HasIn method
205  {
206  // the point is added
207  this->AddSBPoint(*itPt);
208  if(this->GetSize() >= (unsigned int)pMinPts)
209  {
210  return;
211  }
212  // restart from scratch
213  itPt = pPtsToCheck->begin();
214  }else
215  {
216  ++itPt;
217  }
218  }
219 }
220 
221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222 
223 bool ClusterSBPoints::HasIn(const SBPoint* pPtToCheck, G4double pMinDist)
224 {
225  // check if the given point is near one of the cluster's point
226  std::set<SBPoint*>::iterator itClusPt;
227  for(itClusPt = fpRegisteredSBPoints.begin();
228  itClusPt != fpRegisteredSBPoints.end();
229  ++itClusPt)
230  {
231  // if are two different pts
232  if((*pPtToCheck != *(*itClusPt)))
233  {
234  // if close enought of an include point of the cluster
235  if( AreOnTheSameCluster(pPtToCheck, *itClusPt, pMinDist) )
236  {
237  return true;
238  }
239  }
240  }
241  return false;
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245 
247  G4double pMinDist)
248 {
249 
250  G4double x1=pPtToCheck->GetPosition().x()/nm;
251  G4double y1=pPtToCheck->GetPosition().y()/nm;
252  G4double z1=pPtToCheck->GetPosition().z()/nm;
253 
254  G4double x2=this->GetBarycenter().x()/nm;
255  G4double y2=this->GetBarycenter().y()/nm;
256  G4double z2=this->GetBarycenter().z()/nm;
257 
258  // if the two points are closed enough
259  if(((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2))<=
260  (pMinDist/nm*pMinDist/nm))
261  {
262  return true;
263  }else
264  {
265  return false;
266  }
267 }
268 
269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
270 
274 {
275  std::set<SBPoint*> points = pCluster->GetRegistredSBPoints();
276  pCluster->Clear();
277  std::set<SBPoint*>::iterator itPt;
278  for(itPt = points.begin(); itPt != points.end(); ++itPt)
279  {
280  this->AddSBPoint(*itPt);
281  }
282 }
283