source: trunk/examples/advanced/brachytherapy/src/BrachyPhantomSD.cc@ 1211

Last change on this file since 1211 was 807, checked in by garnier, 17 years ago

update

File size: 4.6 KB
Line 
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//
27// Code developed by:
28// S. Agostinelli, F. Foppiano, S. Garelli , M. Tropeano, S.Guatelli
29//
30// ********************************
31// * *
32// * BrachyPhantomSD.cc *
33// * *
34// ********************************
35//
36// $Id: BrachyPhantomSD.cc,v 1.13 2006/06/29 15:48:39 gunter Exp $
37// GEANT4 tag $Name: $
38//
39#include "BrachyPhantomSD.hh"
40#include "BrachyAnalysisManager.hh"
41#include "BrachyDetectorConstruction.hh"
42#include "G4Track.hh"
43#include "G4LogicalVolume.hh"
44#include "G4VPhysicalVolume.hh"
45#include "G4Step.hh"
46#include "G4VTouchable.hh"
47#include "G4TouchableHistory.hh"
48#include "G4SDManager.hh"
49#include "G4ParticleDefinition.hh"
50
51BrachyPhantomSD::BrachyPhantomSD(G4String name):G4VSensitiveDetector(name)
52{
53}
54
55BrachyPhantomSD::~BrachyPhantomSD()
56{
57
58}
59
60void BrachyPhantomSD::Initialize(G4HCofThisEvent*)
61{
62}
63
64G4bool BrachyPhantomSD::ProcessHits(G4Step* aStep, G4TouchableHistory* ROhist)
65{
66 //
67 // The energy deposit of the hits is stored in histograms and ntuples
68 //
69
70 if(!ROhist)
71 return false;
72
73 // Check the volume
74 if(aStep -> GetPreStepPoint() -> GetPhysicalVolume() ->
75 GetName() != "PhantomPhys")
76 return false;
77
78 G4double energyDeposit = aStep -> GetTotalEnergyDeposit();
79
80 // Check that the energy deposit is not null
81 if(energyDeposit == 0.)
82 return false;
83
84 if(energyDeposit != 0)
85 {
86 // Read Voxel indexes:
87 // i is the x index,
88 // j is the y index
89 // k is the z index
90 G4int j = ROhist -> GetReplicaNumber();
91 G4int k = ROhist -> GetReplicaNumber(1);
92 G4int i = ROhist -> GetReplicaNumber(2);
93
94 G4int numberOfVoxel = 300;
95 G4double voxelWidth = 1. *mm;
96
97 // Retrieve the coordinates of the center of the voxel where
98 // the energy deposit is located
99 x = ( - numberOfVoxel + 1+ 2*i )* voxelWidth/2;
100 y = ( - numberOfVoxel + 1+ 2*j )* voxelWidth/2;
101 z = ( - numberOfVoxel + 1+ 2*k )* voxelWidth/2;
102
103#ifdef G4ANALYSIS_USE
104 BrachyAnalysisManager* analysis =
105 BrachyAnalysisManager::getInstance();
106
107 // Fill the ntuple with position and energy deposit in the phantom
108 analysis -> FillNtupleWithEnergy(x,y,z,energyDeposit/MeV);
109
110 if (y < 0.8 * mm && y > -0.8 * mm)
111 {
112 // Fill a 2D histogram with the energy deposit in the plane
113 // containing the source
114 analysis -> FillHistogramWithEnergy(x,z,energyDeposit/MeV);
115
116 // Fill 1D histogram with the energy deposit
117 // along the axis perpendicular to the main axis of the source
118 if (z < 0.8 * mm && z > -0.8 * mm)
119 analysis -> DoseDistribution(x,energyDeposit/MeV);
120 }
121#endif
122 }
123 return true;
124}
125
126void BrachyPhantomSD::EndOfEvent(G4HCofThisEvent*)
127{
128}
129
130void BrachyPhantomSD::clear()
131{
132}
133
134void BrachyPhantomSD::DrawAll()
135{
136}
137
138void BrachyPhantomSD::PrintAll()
139{
140}
141
142
143
Note: See TracBrowser for help on using the repository browser.