Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/exp_microdosimetry/src/SensitiveDetector.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 // ********************************************************************
  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 // Authors: Susanna Guatelli and Francesco Romano
 27 // susanna@uow.edu.au, francesco.romano@ct.infn.it
 28 //
 29 //  code based on the basic example B2
 30 //
 31 #include "SensitiveDetector.hh"
 32 #include "G4HCofThisEvent.hh"
 33 #include "G4Step.hh"
 34 #include "G4ThreeVector.hh"
 35 #include "G4SDManager.hh"
 36 #include "G4ios.hh"
 37 #include "G4SystemOfUnits.hh"
 38 #include "G4RunManager.hh"
 39 
 40 SensitiveDetector::SensitiveDetector(const G4String& name,
 41                          const G4String& hitsCollectionName, 
 42                          G4bool isThisAMicrodosimeter,
 43                          AnalysisManager* analysis_manager) 
 44  : G4VSensitiveDetector(name),
 45    fHitsCollection(NULL)
 46 {
 47   collectionName.insert(hitsCollectionName);
 48   analysis = analysis_manager;
 49   isMicrod = isThisAMicrodosimeter; // only false for the E-stage of the telescope
 50 }
 51 
 52 SensitiveDetector::~SensitiveDetector() 
 53 {}
 54 
 55 void SensitiveDetector::Initialize(G4HCofThisEvent* hce)
 56 {
 57   // Create hits collection
 58 
 59   fHitsCollection 
 60     = new SensitiveDetectorHitsCollection(SensitiveDetectorName, collectionName[0]); 
 61 
 62   // Add this collection in hce
 63 
 64   G4int hcID 
 65     = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
 66   hce->AddHitsCollection( hcID, fHitsCollection ); 
 67 }
 68 
 69 G4bool SensitiveDetector::ProcessHits(G4Step* aStep, 
 70                                      G4TouchableHistory*)
 71 {  
 72   // energy deposit
 73   G4double edep = aStep->GetTotalEnergyDeposit();
 74   
 75   if (edep==0.) return false;
 76 
 77   // Uncomment the following lines to only score specific physical volumes 
 78   /*
 79   G4String volumeName = aStep -> GetPreStepPoint() -> GetPhysicalVolume()-> GetName();
 80   if(volumeName != "SV_phys1" && volumeName != "sen_bridge" && volumeName != "physSensitiveBridgeVolume" ) 
 81     return false;  
 82   */
 83   
 84   SensitiveDetectorHit* newHit = new SensitiveDetectorHit();
 85 
 86   newHit -> SetEdep(edep);
 87   
 88   // step length
 89   G4double len = aStep->GetStepLength();
 90 
 91   // only consider ions, protons, and neutrons for path length
 92   const G4ParticleDefinition* thisParticle = aStep ->
 93                       GetTrack() -> GetParticleDefinition() ;
 94   G4String particleName = thisParticle -> GetParticleName();
 95   G4int particleZ = thisParticle -> GetAtomicNumber();
 96 
 97   if ( (particleZ < 1) && (particleName != "neutron"))
 98     len = 0.;
 99 
100   newHit->SetPath(len);
101   
102   fHitsCollection -> insert( newHit );
103 
104   return true;
105 }
106 
107 void SensitiveDetector::EndOfEvent(G4HCofThisEvent*)
108 {
109 #ifdef ANALYSIS_USE
110 // Initialisation of total energy deposition per event to zero
111  G4double totalEdepInOneEvent=0;
112  G4double totalPathLengthInOneEvent=0;
113  
114  G4int NbHits = fHitsCollection->entries();
115    //G4cout << "number of hits " <<NbHits << G4endl;
116  
117  G4double edep, len; 
118  
119  for (G4int i=0;i<NbHits;i++) 
120  {
121      edep = (*fHitsCollection)[i]->GetEdep();
122      totalEdepInOneEvent = totalEdepInOneEvent + edep;
123  
124      len = (*fHitsCollection)[i]->GetPath();
125      totalPathLengthInOneEvent = totalPathLengthInOneEvent + len;
126  }
127  
128  G4int eventID = G4RunManager::GetRunManager() ->
129      GetCurrentEvent() -> GetEventID();
130 
131 
132  if (totalEdepInOneEvent!=0)
133  {
134      if (isMicrod==true) analysis -> StoreEnergyDeposition(totalEdepInOneEvent/keV, totalPathLengthInOneEvent/um, eventID);
135      else analysis -> StoreSecondStageEnergyDeposition(totalEdepInOneEvent/keV, eventID); 
136  }
137 #endif 
138 }
139