Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ExN05EMShowerModel.cc
Go to the documentation of this file.
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 // $Id$
28 //
29 #include "ExN05EMShowerModel.hh"
30 #include "ExN05EnergySpot.hh"
31 
32 #include "Randomize.hh"
33 
34 #include "G4Electron.hh"
35 #include "G4Positron.hh"
36 #include "G4Gamma.hh"
38 #include "G4VSensitiveDetector.hh"
39 #include "G4TouchableHandle.hh"
40 #include "G4PhysicalConstants.hh"
41 #include "G4SystemOfUnits.hh"
42 
44 : G4VFastSimulationModel(modelName, envelope)
45 {
46  fFakeStep = new G4Step();
47  fFakePreStepPoint = fFakeStep->GetPreStepPoint();
48  fFakePostStepPoint = fFakeStep->GetPostStepPoint();
49  fTouchableHandle = new G4TouchableHistory();
50  fpNavigator = new G4Navigator();
51  fNaviSetup = false;
52 }
53 
55 : G4VFastSimulationModel(modelName)
56 {
57  fFakeStep = new G4Step();
58  fFakePreStepPoint = fFakeStep->GetPreStepPoint();
59  fFakePostStepPoint = fFakeStep->GetPostStepPoint();
60  fTouchableHandle = new G4TouchableHistory();
61  fpNavigator = new G4Navigator();
62  fNaviSetup = false;
63 }
64 
66 {
67  delete fFakeStep;
68  delete fpNavigator;
69 }
70 
72 {
73  return
74  &particleType == G4Electron::ElectronDefinition() ||
75  &particleType == G4Positron::PositronDefinition() ||
76  &particleType == G4Gamma::GammaDefinition();
77 }
78 
80 {
81  // Applies the parameterisation above 100 MeV:
82  return fastTrack.GetPrimaryTrack()->GetKineticEnergy() > 100*MeV;
83 }
84 
85 void ExN05EMShowerModel::DoIt(const G4FastTrack& fastTrack,
86  G4FastStep& fastStep)
87 {
88  // Kill the parameterised particle:
89  fastStep.KillPrimaryTrack();
90  fastStep.ProposePrimaryTrackPathLength(0.0);
92 
93  // split into "energy spots" energy according to the shower shape:
94  Explode(fastTrack);
95 
96  // and put those energy spots into the crystals:
97  BuildDetectorResponse();
98 
99 }
100 
101 void ExN05EMShowerModel::Explode(const G4FastTrack& fastTrack)
102 {
103  //-----------------------------------------------------
104  //
105  //-----------------------------------------------------
106 
107  // Reduced quantities:
108  // -- critical energy in CsI:
109  G4double Ec = 800*MeV/(54. + 1.2); // 54 = mean Z of CsI
110  G4double Energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
111  G4double y = Energy/Ec;
112 
113  // compute value of parameter "a" of longitudinal profile, b assumed = 0.5
114  G4double a, tmax, b(0.5), C;
115  if (fastTrack.GetPrimaryTrack()->GetDefinition() == G4Gamma::GammaDefinition()) C = 0.5;
116  else C = -0.5;
117  tmax = 1.0 * (std::log(y) + C);
118  a = 1.0 + b*tmax;
119 
120  // t : reduced quantity = z/X0:
121  G4double t, bt;
122  G4Material* CsI = G4Material::GetMaterial("CsI");
123  G4double X0 = CsI->GetRadlen();
124  // Moliere radius:
125  G4double Es = 21*MeV;
126  G4double Rm = X0*Es/Ec;
127 
128  // axis of the shower, in global reference frame:
129  G4ThreeVector xShower, yShower, zShower;
130  zShower = fastTrack.GetPrimaryTrack()->GetMomentumDirection();
131  xShower = zShower.orthogonal();
132  yShower = zShower.cross(xShower);
133  // starting point of the shower:
134  G4ThreeVector sShower = fastTrack.GetPrimaryTrack()->GetPosition();
135 
136  // We shoot 100 spots of energy:
137  G4int nSpots = 100;
138  G4double deposit = Energy/double(nSpots);
139  ExN05EnergySpot eSpot;
140  eSpot.SetEnergy(deposit);
141  G4ThreeVector ePoint;
142  G4double z, r, phi;
143 
144  feSpotList.clear();
145  for (int i = 0; i < nSpots; i++)
146  {
147  // Longitudinal profile:
148  // -- shoot z according to Gamma distribution:
149  bt = CLHEP::RandGamma::shoot(a,1.0);
150  t = bt/b;
151  z = t*X0;
152 
153  // transverse profile:
154  // we set 90% of energy in one Rm,
155  // the rest between 1 and 3.5 Rm:
156  G4double xr = G4UniformRand();
157  if (xr < 0.9) r = xr/0.9*Rm;
158  else r = ((xr - 0.9)/0.1*2.5 + 1.0)*Rm;
159  phi = G4UniformRand()*twopi;
160 
161  // build the position:
162  ePoint = sShower +
163  z*zShower +
164  r*std::cos(phi)*xShower + r*std::sin(phi)*yShower;
165 
166  // and the energy spot:
167  eSpot.SetPosition(ePoint);
168 
169  // Records the eSpot:
170  feSpotList.push_back(eSpot);
171  }
172 }
173 
174 
175 void ExN05EMShowerModel::BuildDetectorResponse()
176 {
177  // Does the assignation of the energy spots to the sensitive volumes:
178  for (size_t i = 0; i < feSpotList.size(); i++)
179  {
180  // Draw the energy spot:
181  feSpotList[i].Draw();
182  // feSpotList[i].Print();
183 
184  // "converts" the energy spot into the fake
185  // G4Step to pass to sensitive detector:
186  AssignSpotAndCallHit(feSpotList[i]);
187  }
188 }
189 
190 
191 void ExN05EMShowerModel::AssignSpotAndCallHit(const ExN05EnergySpot &eSpot)
192 {
193  //
194  // "converts" the energy spot into the fake
195  // G4Step to pass to sensitive detector:
196  //
197  FillFakeStep(eSpot);
198 
199  //
200  // call sensitive part: taken/adapted from the stepping:
201  // Send G4Step information to Hit/Dig if the volume is sensitive
202  //
203  G4VPhysicalVolume* pCurrentVolume =
204  fFakeStep->GetPreStepPoint()->GetPhysicalVolume();
205  G4VSensitiveDetector* pSensitive;
206 
207  if( pCurrentVolume != 0 )
208  {
209  pSensitive = pCurrentVolume->GetLogicalVolume()->
210  GetSensitiveDetector();
211  if( pSensitive != 0 )
212  {
213  pSensitive->Hit(fFakeStep);
214  }
215  }
216 }
217 
218 
219 void ExN05EMShowerModel::FillFakeStep(const ExN05EnergySpot &eSpot)
220 {
221  //-----------------------------------------------------------
222  // find in which volume the spot is.
223  //-----------------------------------------------------------
224  if (!fNaviSetup)
225  {
226  fpNavigator->
228  GetNavigatorForTracking()->GetWorldVolume());
229  fpNavigator->
230  LocateGlobalPointAndUpdateTouchableHandle(eSpot.GetPosition(),
231  G4ThreeVector(0.,0.,0.),
232  fTouchableHandle,
233  false);
234  fNaviSetup = true;
235  }
236  else
237  {
238  fpNavigator->
239  LocateGlobalPointAndUpdateTouchableHandle(eSpot.GetPosition(),
240  G4ThreeVector(0.,0.,0.),
241  fTouchableHandle);
242  }
243  //--------------------------------------
244  // Fills attribute of the G4Step needed
245  // by our sensitive detector:
246  //-------------------------------------
247  // set touchable volume at PreStepPoint:
248  fFakePreStepPoint->SetTouchableHandle(fTouchableHandle);
249  // set total energy deposit:
250  fFakeStep->SetTotalEnergyDeposit(eSpot.GetEnergy());
251 }