Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Data Fields
HadrontherapyDetectorSD Class Reference

#include <HadrontherapyDetectorSD.hh>

Inheritance diagram for HadrontherapyDetectorSD:
G4VSensitiveDetector

Public Member Functions

 HadrontherapyDetectorSD (G4String name)
 
 ~HadrontherapyDetectorSD ()
 
void Initialize (G4HCofThisEvent *)
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist)
 
void EndOfEvent (G4HCofThisEvent *HCE)
 
- Public Member Functions inherited from G4VSensitiveDetector
 G4VSensitiveDetector (G4String name)
 
 G4VSensitiveDetector (const G4VSensitiveDetector &right)
 
virtual ~G4VSensitiveDetector ()
 
const 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
 

Data Fields

std::ofstream ofs
 

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

Definition at line 45 of file HadrontherapyDetectorSD.hh.

Constructor & Destructor Documentation

HadrontherapyDetectorSD::HadrontherapyDetectorSD ( G4String  name)

Definition at line 49 of file HadrontherapyDetectorSD.cc.

References G4VSensitiveDetector::collectionName, and G4CollectionNameVector::insert().

49  :
51 {
52  G4String HCname;
53  collectionName.insert(HCname="HadrontherapyDetectorHitsCollection");
54  HitsCollection = NULL;
55  sensitiveDetectorName = name;
56 
57 }
const XML_Char * name
G4VSensitiveDetector(G4String name)
G4CollectionNameVector collectionName
HadrontherapyDetectorSD::~HadrontherapyDetectorSD ( )

Definition at line 60 of file HadrontherapyDetectorSD.cc.

61 {
62 }

Member Function Documentation

void HadrontherapyDetectorSD::EndOfEvent ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VSensitiveDetector.

Definition at line 232 of file HadrontherapyDetectorSD.cc.

References G4VSensitiveDetector::GetCollectionID().

233 {
234 
235  static G4int HCID = -1;
236  if(HCID < 0)
237  {
238  HCID = GetCollectionID(0);
239  }
240 
241  HCE -> AddHitsCollection(HCID,HitsCollection);
242 }
int G4int
Definition: G4Types.hh:78
virtual G4int GetCollectionID(G4int i)
void HadrontherapyDetectorSD::Initialize ( G4HCofThisEvent )
virtual

Reimplemented from G4VSensitiveDetector.

Definition at line 65 of file HadrontherapyDetectorSD.cc.

References G4VSensitiveDetector::collectionName.

66 {
67 
68  HitsCollection = new HadrontherapyDetectorHitsCollection(sensitiveDetectorName,
69  collectionName[0]);
70 }
G4THitsCollection< HadrontherapyDetectorHit > HadrontherapyDetectorHitsCollection
G4CollectionNameVector collectionName
G4bool HadrontherapyDetectorSD::ProcessHits ( G4Step aStep,
G4TouchableHistory ROhist 
)
virtual

Implements G4VSensitiveDetector.

Definition at line 73 of file HadrontherapyDetectorSD.cc.

References HadrontherapyLet::GetInstance(), HadrontherapyAnalysisManager::GetInstance(), HadrontherapyMatrix::GetInstance(), G4VSensitiveDetector::GetName(), G4ParticleDefinition::GetPDGEncoding(), G4Step::GetPreStepPoint(), G4VTouchable::GetReplicaNumber(), G4StepPoint::GetTouchable(), python.hepunit::keV, and python.hepunit::MeV.

74 {
75 
76 
77  if (aStep -> GetPreStepPoint() -> GetPhysicalVolume() -> GetName() != "RODetectorZDivisionPhys") return false;
78 
79 
80  // Get kinetic energy
81  G4Track * theTrack = aStep -> GetTrack();
82  G4double kineticEnergy = theTrack -> GetKineticEnergy();
83 
84  G4ParticleDefinition *particleDef = theTrack -> GetDefinition();
85  //Get particle name
86  G4String particleName = particleDef -> GetParticleName();
87 
88  G4int pdg = particleDef ->GetPDGEncoding();
89 
90  // Get unique track_id (in an event)
91  G4int trackID = theTrack -> GetTrackID();
92 
93  G4double energyDeposit = aStep -> GetTotalEnergyDeposit();
94 
95  G4double DX = aStep -> GetStepLength();
96  G4int Z = particleDef-> GetAtomicNumber();
97  G4int A = particleDef-> GetAtomicMass();
98 
99 
100  // Read voxel indexes: i is the x index, k is the z index
101  const G4VTouchable* touchable = aStep->GetPreStepPoint()->GetTouchable();
102  G4int k = touchable->GetReplicaNumber(0);
103  G4int i = touchable->GetReplicaNumber(2);
104  G4int j = touchable->GetReplicaNumber(1);
105 
106 
107 
108 #ifdef G4ANALYSIS_USE_ROOT
110 #endif
111 
114 
115 
116  // ******************** let ***************************
117  if (let)
118  {
119  if ( !(Z==0 && A==1) ) // All but not neutrons
120  {
121  if( energyDeposit>0. && DX >0. )
122  {
123  if (pdg !=22) // not gamma
124  {
125  let -> FillEnergySpectrum(trackID, particleDef,energyDeposit, DX, i, j, k);
126  }
127  else if (kineticEnergy > 50.*keV) // gamma cut
128  {
129  let -> FillEnergySpectrum(trackID, particleDef,energyDeposit, DX, i , j, k);
130  }
131 
132  }
133  }
134  }
135 
136 
137 
138  // ******************** let ***************************
139 
140 
141  if (matrix)
142  {
143 
144  // Increment Fluences & accumulate energy spectra
145  // Hit voxels are marked with track_id throught hitTrack matrix
146  G4int* hitTrack = matrix -> GetHitTrack(i,j,k); // hitTrack MUST BE cleared at every eventAction!
147  if ( *hitTrack != trackID )
148  {
149  *hitTrack = trackID;
150  /*
151  * Fill FLUENCE data for every single nuclide
152  * Exclude e-, neutrons, gamma, ...
153  */
154  if ( Z >= 1) matrix -> Fill(trackID, particleDef, i, j, k, 0, true);
155 
156 
157 #ifdef G4ANALYSIS_USE_ROOT
158 /*
159  // Fragments kinetic energy (ntuple)
160  if (trackID !=1 && Z>=1)
161  {
162  // First step kinetic energy for every fragment
163  analysis -> FillKineticFragmentTuple(i, j, k, A, Z, kineticEnergy/MeV);
164  }
165  // Kinetic energy spectra for primary particles
166 
167  if ( trackID == 1 && i == 0)
168  {
169  // First step kinetic energy for primaries only
170  analysis -> FillKineticEnergyPrimaryNTuple(i, j, k, kineticEnergy/MeV);
171  }
172 */
173 #endif
174  }
175 
176  if(energyDeposit != 0)
177  {
178 /*
179  * This method will fill a dose matrix for every single nuclide.
180  * A method of the HadrontherapyMatrix class (StoreDoseFluenceAscii())
181  * is called automatically at the end of main (or via the macro command /analysis/writeDoseFile.
182  * It permits to store all dose/fluence data into a single plane ASCII file.
183 */
184  // if (A==1 && Z==1) // primary and sec. protons
185  if ( Z>=1 ) // exclude e-, neutrons, gamma, ...
186  matrix -> Fill(trackID, particleDef, i, j, k, energyDeposit);
187  /*
188  * Create a hit with the information of position is in the detector
189  */
191  detectorHit -> SetEdepAndPosition(i, j, k, energyDeposit);
192  HitsCollection -> insert(detectorHit);
193  }
194  }
195 
196 #ifdef G4ANALYSIS_USE_ROOT
197  if(energyDeposit != 0)
198  {
199  if(trackID != 1)
200  {
201  if (particleName == "proton")
202  analysis -> SecondaryProtonEnergyDeposit(i, energyDeposit/MeV);
203 
204  else if (particleName == "neutron")
205  analysis -> SecondaryNeutronEnergyDeposit(i, energyDeposit/MeV);
206 
207  else if (particleName == "alpha")
208  analysis -> SecondaryAlphaEnergyDeposit(i, energyDeposit/MeV);
209 
210  else if (particleName == "gamma")
211  analysis -> SecondaryGammaEnergyDeposit(i, energyDeposit/MeV);
212 
213  else if (particleName == "e-")
214  analysis -> SecondaryElectronEnergyDeposit(i, energyDeposit/MeV);
215 
216  else if (particleName == "triton")
217  analysis -> SecondaryTritonEnergyDeposit(i, energyDeposit/MeV);
218 
219  else if (particleName == "deuteron")
220  analysis -> SecondaryDeuteronEnergyDeposit(i, energyDeposit/MeV);
221 
222  else if (particleName == "pi+" || particleName == "pi-" || particleName == "pi0")
223  analysis -> SecondaryPionEnergyDeposit(i, energyDeposit/MeV);
224  }
225  }
226 #endif
227 
228  return true;
229 }
static HadrontherapyAnalysisManager * GetInstance()
static HadrontherapyLet * GetInstance()
const G4VTouchable * GetTouchable() const
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
static HadrontherapyMatrix * GetInstance()
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
double G4double
Definition: G4Types.hh:76

Field Documentation

std::ofstream HadrontherapyDetectorSD::ofs

Definition at line 51 of file HadrontherapyDetectorSD.hh.


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