Geant4  10.02.p03
TrackerSD Class Reference

#include <TrackerSD.hh>

Inheritance diagram for TrackerSD:
Collaboration diagram for TrackerSD:

Public Member Functions

 TrackerSD (const G4String &name, const G4String &hitsCollectionName)
 
virtual ~TrackerSD ()
 
virtual void Initialize (G4HCofThisEvent *hitCollection)
 
virtual G4bool ProcessHits (G4Step *step, G4TouchableHistory *history)
 
virtual void EndOfEvent (G4HCofThisEvent *hitCollection)
 
- Public Member Functions inherited from G4VSensitiveDetector
 G4VSensitiveDetector (G4String name)
 
 G4VSensitiveDetector (const G4VSensitiveDetector &right)
 
virtual ~G4VSensitiveDetector ()
 
G4VSensitiveDetectoroperator= (const G4VSensitiveDetector &right)
 
G4int operator== (const G4VSensitiveDetector &right) const
 
G4int operator!= (const G4VSensitiveDetector &right) const
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
G4bool Hit (G4Step *aStep)
 
void SetROgeometry (G4VReadOutGeometry *value)
 
void SetFilter (G4VSDFilter *value)
 
G4int GetNumberOfCollections () const
 
G4String GetCollectionName (G4int id) const
 
void SetVerboseLevel (G4int vl)
 
void Activate (G4bool activeFlag)
 
G4bool isActive () const
 
G4String GetName () const
 
G4String GetPathName () const
 
G4String GetFullPathName () const
 
G4VReadOutGeometryGetROgeometry () const
 
G4VSDFilterGetFilter () const
 
virtual G4VSensitiveDetectorClone () const
 

Private Attributes

TrackerHitsCollectionfHitsCollection
 

Additional Inherited Members

- Protected Member Functions inherited from G4VSensitiveDetector
virtual G4int GetCollectionID (G4int i)
 
- Protected Attributes inherited from G4VSensitiveDetector
G4CollectionNameVector collectionName
 
G4String SensitiveDetectorName
 
G4String thePathName
 
G4String fullPathName
 
G4int verboseLevel
 
G4bool active
 
G4VReadOutGeometryROgeometry
 
G4VSDFilterfilter
 

Detailed Description

Tracker sensitive detector class

The hits are accounted in hits in ProcessHits() function which is called by Geant4 kernel at each step. A hit is created with each step with non zero energy deposit.

Definition at line 52 of file TrackerSD.hh.

Constructor & Destructor Documentation

◆ TrackerSD()

TrackerSD::TrackerSD ( const G4String name,
const G4String hitsCollectionName 
)

Definition at line 47 of file TrackerSD.cc.

49  : G4VSensitiveDetector(name),
50  fHitsCollection(NULL)
51 {
52  collectionName.insert(hitsCollectionName);
53 }
TrackerHitsCollection * fHitsCollection
Definition: TrackerSD.hh:65
G4VSensitiveDetector(G4String name)
G4CollectionNameVector collectionName
Here is the call graph for this function:

◆ ~TrackerSD()

TrackerSD::~TrackerSD ( )
virtual

Definition at line 57 of file TrackerSD.cc.

58 {}

Member Function Documentation

◆ EndOfEvent()

void TrackerSD::EndOfEvent ( G4HCofThisEvent hitCollection)
virtual

Reimplemented from G4VSensitiveDetector.

Definition at line 101 of file TrackerSD.cc.

102 {
103  G4int nofHits = fHitsCollection->entries();
104 
105  if ( verboseLevel>1 )
106 
107  G4cout << G4endl
108  << "-------->Hits Collection: in this event they are "
109  << nofHits
110  << " hits in the target volume " << G4endl;
111 
112  // PROCESSING OF MICRODOSIMETRY Y & Z SPECTRA
113 
114  // main loop
115 
116  // *************************************
117  // Please select herebelow :
118  // 1 - the number of sampling per incident envent
119  // variable name = nbSampling
120  // (it is set to 100 by default)
121  //
122  //
123  // 2 - the radius of the target sphere
124  // variable name = radius
125  // (it is set to 50 nm by default)
126  //
127 
128  G4int nbSampling = 1000;
129  G4double radius = 1.50*um;
130 
131  // ************************************
132 
133  for ( G4int k=0; k<nbSampling; k++ )
134  {
135 
136  //***************
137  // FIRST SAMPLING: Y
138  //***************
139 
140  // select random hit
141  G4int randHit=0; // Runs from 0 to number of hits - 1
142  randHit = static_cast<G4int>( G4UniformRand()*nofHits );
143 
144  //G4cout << static_cast<G4int>(2.7)
145  //<< "======> random selection of hit number randHit ="
146  // << randHit << G4endl;
147 
148  // get random hit position
149  G4ThreeVector hitPos = (*fHitsCollection)[randHit]->GetPos();
150  //G4cout << "======> random hit position x/nm =" << hitPos.x()/nm << G4endl;
151  //G4cout << "======> random hit position y/nm =" << hitPos.y()/nm << G4endl;
152  //G4cout << "======> random hit position z/nm =" << hitPos.z()/nm << G4endl;
153 
154  // select random center of sphere within radius
155  G4double chord = 4.*radius/3;
156  //G4double density = 1 * g/cm3;
157  //G4double mass = (4./3)*CLHEP::pi*radius*radius*radius*density;
158 
159  G4ThreeVector randDir = G4RandomDirection();
160  G4double randRadius = G4UniformRand()*radius;
161  G4ThreeVector randCenterPos = randRadius*randDir + hitPos;
162 
163  // search for neighbouring hits in the sphere and cumulate deposited energy
164  // in epsilon
165  G4double epsilon = 0;
166  G4int nbEdep = 0;
167 
168  for ( G4int i=0; i<nofHits; i++ )
169  {
170  G4ThreeVector localPos = (*fHitsCollection)[i]->GetPos();;
171  if (
172  (localPos.x()-randCenterPos.x()) * (localPos.x()-randCenterPos.x()) +
173  (localPos.y()-randCenterPos.y()) * (localPos.y()-randCenterPos.y()) +
174  (localPos.z()-randCenterPos.z()) * (localPos.z()-randCenterPos.z())
175  <= radius*radius
176  )
177 
178  {
179  epsilon = epsilon + (*fHitsCollection)[i]->GetEdep() ;
180  nbEdep = nbEdep+1;
181  }
182 
183  }
184 
185 
186  //***************
187 
188  /*
189  // for testing only
190 
191  G4cout << "======> for hit number #" << randHit <<
192  ", we collect "
193  << nbEdep << " energy depositions in a sphere of radius "
194  << radius/nm << " nm and mass "
195  << mass/kg << " kg for a total of "
196  << epsilon/eV << " eV or "
197  << (epsilon/joule)/(mass/kg) << " Gy" << G4endl;
198  G4cout << "-" << G4endl;
199 
200  FILE* myFile;
201  myFile=fopen("yz.txt","a");
202  fprintf(myFile,"%e %e %e\n",radius/nm,(epsilon/eV)/(chord/nm),
203  (epsilon/joule)/(mass/kg));
204  fclose(myFile);
205  */
206 
207  // get analysis manager
208  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
209 
210  // fill ntuple
211  analysisManager->FillNtupleDColumn(0, radius/nm);
212  analysisManager->FillNtupleDColumn(1, (epsilon/eV)/(chord/nm));
213  analysisManager->AddNtupleRow();
214 
215  } // END MAIN LOOP ON SAMPLING
216 
217 }
G4ThreeVector G4RandomDirection()
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static const double nm
Definition: G4SIunits.hh:111
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:66
double x() const
TrackerHitsCollection * fHitsCollection
Definition: TrackerSD.hh:65
static const double eV
Definition: G4SIunits.hh:212
double y() const
double z() const
static const double um
Definition: G4SIunits.hh:112
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
double epsilon(double density, double temperature)
Here is the call graph for this function:

◆ Initialize()

void TrackerSD::Initialize ( G4HCofThisEvent hitCollection)
virtual

Reimplemented from G4VSensitiveDetector.

Definition at line 62 of file TrackerSD.cc.

63 {
64  // Create hits collection
65 
68 
69  // Add this collection in hce
70 
71  G4int hcID
73  hce->AddHitsCollection( hcID, fHitsCollection );
74 }
G4THitsCollection< TrackerHit > TrackerHitsCollection
Definition: TrackerHit.hh:87
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
int G4int
Definition: G4Types.hh:78
TrackerHitsCollection * fHitsCollection
Definition: TrackerSD.hh:65
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
G4CollectionNameVector collectionName
Here is the call graph for this function:

◆ ProcessHits()

G4bool TrackerSD::ProcessHits ( G4Step *  step,
G4TouchableHistory history 
)
virtual

Implements G4VSensitiveDetector.

Definition at line 78 of file TrackerSD.cc.

80 {
81  // energy deposit
82  G4double edep = aStep->GetTotalEnergyDeposit();
83 
84  if (edep==0.) return false;
85 
86  TrackerHit* newHit = new TrackerHit();
87 
88  newHit->SetTrackID (aStep->GetTrack()->GetTrackID());
89  newHit->SetEdep(edep);
90  newHit->SetPos (aStep->GetPostStepPoint()->GetPosition());
91 
92  fHitsCollection->insert( newHit );
93 
94  //newHit->Print();
95 
96  return true;
97 }
void SetTrackID(G4int track)
Definition: TrackerHit.hh:69
Double_t edep
void SetPos(G4ThreeVector xyz)
Definition: TrackerHit.hh:71
TrackerHitsCollection * fHitsCollection
Definition: TrackerSD.hh:65
void SetEdep(G4double de)
Definition: TrackerHit.hh:70
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

Member Data Documentation

◆ fHitsCollection

TrackerHitsCollection* TrackerSD::fHitsCollection
private

Definition at line 65 of file TrackerSD.hh.


The documentation for this class was generated from the following files: