Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:21:51

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // This example is provided by the Geant4-DNA collaboration
0027 // Any report or published results obtained using the Geant4-DNA software
0028 // shall cite the following Geant4-DNA collaboration publication:
0029 // Med. Phys. 37 (2010) 4692-4708
0030 // The Geant4-DNA web site is available at http://geant4-dna.org
0031 //
0032 // Authors: Henri Payno and Yann Perrot
0033 //
0034 //
0035 /// \file ClusteringAlgo.cc
0036 /// \brief Implementation of the ClustreringAlgo class
0037 
0038 #include "ClusteringAlgo.hh"
0039 
0040 #include "ClusteringAlgoMessenger.hh"
0041 
0042 #include "G4SystemOfUnits.hh"
0043 #include "Randomize.hh"
0044 
0045 #include <map>
0046 
0047 using namespace std;
0048 
0049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0050 
0051 ClusteringAlgo::ClusteringAlgo(G4double pEps, G4int pMinPts, G4double pSPointsProb,
0052                                G4double pEMinDamage, G4double pEMaxDamage)
0053   : fEps(pEps),
0054     fMinPts(pMinPts),
0055     fSPointsProb(pSPointsProb),
0056     fEMinDamage(pEMinDamage),
0057     fEMaxDamage(pEMaxDamage)
0058 {
0059   fNextSBPointID = 0;
0060   fpClustAlgoMessenger = new ClusteringAlgoMessenger(this);
0061 }
0062 
0063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0064 
0065 ClusteringAlgo::~ClusteringAlgo()
0066 {
0067   delete fpClustAlgoMessenger;
0068   Purge();
0069 }
0070 
0071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0072 
0073 // Random sampling in space
0074 G4bool ClusteringAlgo::IsInSensitiveArea()
0075 {
0076   return fSPointsProb > G4UniformRand();
0077 }
0078 
0079 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0080 
0081 // Random sampling in energy
0082 G4bool ClusteringAlgo::IsEdepSufficient(G4double pEdep)
0083 {
0084   if (pEdep < fEMinDamage) {
0085     return false;
0086   }
0087 
0088   else if (pEdep > fEMaxDamage) {
0089     return true;
0090   }
0091   else {
0092     G4double proba = (pEdep / eV - fEMinDamage / eV) / (fEMaxDamage / eV - fEMinDamage / eV);
0093     return (proba > G4UniformRand());
0094   }
0095 }
0096 
0097 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0098 
0099 // Add an event interaction to the unregistered damage if
0100 // good conditions (pos and energy) are met
0101 //
0102 
0103 void ClusteringAlgo::RegisterDamage(G4ThreeVector pPos, G4double pEdep)
0104 {
0105   if (IsEdepSufficient(pEdep)) {
0106     if (IsInSensitiveArea()) {
0107       fpSetOfPoints.push_back(new SBPoint(fNextSBPointID++, pPos, pEdep));
0108     }
0109   }
0110 }
0111 
0112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0113 
0114 map<G4int, G4int> ClusteringAlgo::RunClustering()
0115 {
0116   // quick sort style
0117   // create cluster
0118   std::vector<SBPoint*>::iterator itVisitorPt, itObservedPt;
0119   for (itVisitorPt = fpSetOfPoints.begin(); itVisitorPt != fpSetOfPoints.end(); ++itVisitorPt) {
0120     itObservedPt = itVisitorPt;
0121     itObservedPt++;
0122     while (itObservedPt != fpSetOfPoints.end()) {
0123       // if at least one of the two points has not a cluster
0124       if (!((*itObservedPt)->HasCluster() && (*itVisitorPt)->HasCluster())) {
0125         if (AreOnTheSameCluster((*itObservedPt)->GetPosition(), (*itVisitorPt)->GetPosition(),
0126                                 fEps))
0127         {
0128           // if none has a cluster. Create a new one
0129           if (!(*itObservedPt)->HasCluster() && !(*itVisitorPt)->HasCluster()) {
0130             // create the new cluster
0131             set<SBPoint*> clusterPoints;
0132             clusterPoints.insert((*itObservedPt));
0133             clusterPoints.insert((*itVisitorPt));
0134             ClusterSBPoints* lCluster = new ClusterSBPoints(clusterPoints);
0135             assert(lCluster);
0136             fpClusters.push_back(lCluster);
0137             assert(lCluster);
0138             // inform SB point that they are part of a cluster now
0139             assert(lCluster);
0140             (*itObservedPt)->SetCluster(lCluster);
0141             assert(lCluster);
0142             (*itVisitorPt)->SetCluster(lCluster);
0143           }
0144           else {
0145             // add the point to the existing cluster
0146             if ((*itObservedPt)->HasCluster()) {
0147               (*itObservedPt)->GetCluster()->AddSBPoint((*itVisitorPt));
0148               (*itVisitorPt)->SetCluster((*itObservedPt)->GetCluster());
0149             }
0150 
0151             if ((*itVisitorPt)->HasCluster()) {
0152               (*itVisitorPt)->GetCluster()->AddSBPoint((*itObservedPt));
0153               (*itObservedPt)->SetCluster((*itVisitorPt)->GetCluster());
0154             }
0155           }
0156         }
0157       }
0158       ++itObservedPt;
0159     }
0160   }
0161 
0162   // associate isolated points and merge clusters
0163   IncludeUnassociatedPoints();
0164   MergeClusters();
0165 
0166   // return cluster size distribution
0167   return GetClusterSizeDistribution();
0168 }
0169 
0170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0171 
0172 // try to merge cluster between them, based on the distance between barycenters
0173 void ClusteringAlgo::MergeClusters()
0174 {
0175   std::vector<ClusterSBPoints*>::iterator itCluster1, itCluster2;
0176   for (itCluster1 = fpClusters.begin(); itCluster1 != fpClusters.end(); ++itCluster1) {
0177     G4ThreeVector baryCenterClust1 = (*itCluster1)->GetBarycenter();
0178     itCluster2 = itCluster1;
0179     itCluster2++;
0180     while (itCluster2 != fpClusters.end()) {
0181       G4ThreeVector baryCenterClust2 = (*itCluster2)->GetBarycenter();
0182       // if we can merge both cluster
0183       if (AreOnTheSameCluster(baryCenterClust1, baryCenterClust2, fEps)) {
0184         (*itCluster1)->MergeWith(*itCluster2);
0185         delete *itCluster2;
0186         fpClusters.erase(itCluster2);
0187         return MergeClusters();
0188       }
0189       else {
0190         itCluster2++;
0191       }
0192     }
0193   }
0194 }
0195 
0196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0197 
0198 void ClusteringAlgo::IncludeUnassociatedPoints()
0199 {
0200   std::vector<SBPoint*>::iterator itVisitorPt;
0201   // Associate all point not in a cluster if possible ( to the first found cluster)
0202   for (itVisitorPt = fpSetOfPoints.begin(); itVisitorPt != fpSetOfPoints.end(); ++itVisitorPt) {
0203     if (!(*itVisitorPt)->HasCluster()) {
0204       FindCluster(*itVisitorPt);
0205     }
0206   }
0207 }
0208 
0209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0210 
0211 bool ClusteringAlgo::FindCluster(SBPoint* pPt)
0212 {
0213   assert(!pPt->HasCluster());
0214   std::vector<ClusterSBPoints*>::iterator itCluster;
0215   for (itCluster = fpClusters.begin(); itCluster != fpClusters.end(); ++itCluster) {
0216     // if((*itCluster)->hasIn(pPt, fEps))
0217     if ((*itCluster)->HasInBarycenter(pPt, fEps)) {
0218       (*itCluster)->AddSBPoint(pPt);
0219       return true;
0220     }
0221   }
0222   return false;
0223 }
0224 
0225 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0226 
0227 bool ClusteringAlgo::AreOnTheSameCluster(G4ThreeVector pPt1, G4ThreeVector pPt2, G4double pMinDist)
0228 {
0229   G4double x1 = pPt1.x() / nm;
0230   G4double y1 = pPt1.y() / nm;
0231   G4double z1 = pPt1.z() / nm;
0232 
0233   G4double x2 = pPt2.x() / nm;
0234   G4double y2 = pPt2.y() / nm;
0235   G4double z2 = pPt2.z() / nm;
0236 
0237   // if the two points are closed enough
0238   if (((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2) + (z1 - z2) * (z1 - z2))
0239       <= (pMinDist / nm * pMinDist / nm))
0240   {
0241     return true;
0242   }
0243   else {
0244     return false;
0245   }
0246 }
0247 
0248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0249 
0250 G4int ClusteringAlgo::GetSSB() const
0251 {
0252   G4int nbSSB = 0;
0253   std::vector<SBPoint*>::const_iterator itSDSPt;
0254   for (itSDSPt = fpSetOfPoints.begin(); itSDSPt != fpSetOfPoints.end(); ++itSDSPt) {
0255     if (!(*itSDSPt)->HasCluster()) {
0256       nbSSB++;
0257     }
0258   }
0259   return nbSSB;
0260 }
0261 
0262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0263 
0264 G4int ClusteringAlgo::GetComplexSSB() const
0265 {
0266   G4int nbSSB = 0;
0267   std::vector<ClusterSBPoints*>::const_iterator itCluster;
0268   for (itCluster = fpClusters.begin(); itCluster != fpClusters.end(); ++itCluster) {
0269     if ((*itCluster)->IsSSB()) {
0270       nbSSB++;
0271     }
0272   }
0273   return nbSSB;
0274 }
0275 
0276 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0277 
0278 G4int ClusteringAlgo::GetDSB() const
0279 {
0280   G4int nbDSB = 0;
0281   std::vector<ClusterSBPoints*>::const_iterator itCluster;
0282   for (itCluster = fpClusters.begin(); itCluster != fpClusters.end(); ++itCluster) {
0283     if ((*itCluster)->IsDSB()) {
0284       nbDSB++;
0285     }
0286   }
0287   return nbDSB;
0288 }
0289 
0290 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0291 
0292 map<G4int, G4int> ClusteringAlgo::GetClusterSizeDistribution()
0293 {
0294   std::map<G4int, G4int> sizeDistribution;
0295   sizeDistribution[1] = GetSSB();
0296   std::vector<ClusterSBPoints*>::const_iterator itCluster;
0297   for (itCluster = fpClusters.begin(); itCluster != fpClusters.end(); itCluster++) {
0298     sizeDistribution[(*itCluster)->GetSize()]++;
0299   }
0300   return sizeDistribution;
0301 }
0302 
0303 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0304 
0305 void ClusteringAlgo::Purge()
0306 {
0307   fNextSBPointID = 0;
0308   std::vector<ClusterSBPoints*>::iterator itCluster;
0309   for (itCluster = fpClusters.begin(); itCluster != fpClusters.end(); ++itCluster) {
0310     delete *itCluster;
0311     *itCluster = NULL;
0312   }
0313   fpClusters.clear();
0314   std::vector<SBPoint*>::iterator itPt;
0315   for (itPt = fpSetOfPoints.begin(); itPt != fpSetOfPoints.end(); ++itPt) {
0316     delete *itPt;
0317     *itPt = NULL;
0318   }
0319   fpSetOfPoints.clear();
0320 }