test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LaserSteppingAction.cc
Go to the documentation of this file.
1 
10 #include "G4ParticleTypes.hh"
11 
13 
15  : theDebugLevel(theConf.getUntrackedParameter<int>("DebugLevel",0)),
16  theEnergyLossScalingFactor(theConf.getUntrackedParameter<double>("EnergyLossScalingFactor",1.0))
17 {
18 }
19 
21 {
22 }
23 
24 void LaserSteppingAction::UserSteppingAction(const G4Step * myStep)
25 {
26  G4Step * theStep = const_cast<G4Step*>(myStep);
27 
28  G4Track * theTrack = theStep->GetTrack();
29 
30  // some debug info
31  {
32  G4TrackStatus isGood = theTrack->GetTrackStatus();
33 
34  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: AC1CMS: the PreStep Material = "
35  << theStep->GetPreStepPoint()->GetMaterial()->GetName()
36  << "\n<LaserSteppingAction::UserSteppingAction(const G4Step *)>: AC1CMS: The Track Status = " << isGood;
37  if ( isGood == fStopAndKill )
38  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: AC1CMS: Track Status = fStopAndKill ";
39 
40  if ( theStep->GetPreStepPoint()->GetProcessDefinedStep() != 0 )
41  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: AC1CMS: PreStep Process = "
42  << theStep->GetPreStepPoint()->GetProcessDefinedStep()->GetProcessName();
43  if ( theStep->GetPostStepPoint()->GetProcessDefinedStep() != 0 )
44  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: AC1CMS: PostStep Process = "
45  << theStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
46  }
47 
48  // ***********************************************************************************************************
49  // Set the EnergyDeposit if the photon is absorbed by a active sensor
50  if ( ( theStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName()== "OpAbsorption" ) )
51  {
52  LogDebug("LaserAlignmentStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step*)>: Photon was absorbed! ";
53 
54  if ( theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector() )
55  {
56  double EnergyLoss = theEnergyLossScalingFactor * theTrack->GetTotalEnergy();
57 
58  // use different energy deposit for the discs depending on the z-position to simulate the variable laser power
59  // Disc 1 TEC2TEC
60  if ( ( (theStep->GetPreStepPoint()->GetPosition().z() > 1262.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1382.5 )
61  || (theStep->GetPreStepPoint()->GetPosition().z() < -1262.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1382.5 ) )
62  && ( ( ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.285 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.295 )
63  || ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.84 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.86 )
64  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 3.63 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 3.66 )
65  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.20 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.23 )
66  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.76 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.80 ) ) ) )
67  { theStep->AddTotalEnergyDeposit(EnergyLoss); } // Disc 1 TEC2TEC
68  // Disc 1
69  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 1262.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1382.5 )
70  || (theStep->GetPreStepPoint()->GetPosition().z() < -1262.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1382.5 ) )
71  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2*0.2*0.2)); } // Disc 1
72  // Disc 2 TEC2TEC
73  else if ( ( (theStep->GetPreStepPoint()->GetPosition().z() > 1402.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1522.5 )
74  || (theStep->GetPreStepPoint()->GetPosition().z() < -1402.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1522.5 ) )
75  && ( ( ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.285 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.295 )
76  || ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.84 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.86 )
77  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 3.63 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 3.66 )
78  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.20 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.23 )
79  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.76 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.80 ) ) ) )
80  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2)); } // Disc 2 TEC2TEC
81  // Disc 2
82  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 1402.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1522.5 )
83  || (theStep->GetPreStepPoint()->GetPosition().z() < -1402.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1522.5 ) )
84  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2*0.2)); } // Disc 2
85  // Disc 3 TEC2TEC
86  else if ( ( (theStep->GetPreStepPoint()->GetPosition().z() > 1542.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1662.5 )
87  || (theStep->GetPreStepPoint()->GetPosition().z() < -1542.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1662.5 ) )
88  && ( ( ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.285 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.295 )
89  || ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.84 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.86 )
90  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 3.63 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 3.66 )
91  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.20 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.23 )
92  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.76 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.80 ) ) ) )
93  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2)); } // Disc 3 TEC2TEC
94  // Disc 3
95  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 1542.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1662.5 )
96  || (theStep->GetPreStepPoint()->GetPosition().z() < -1542.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1662.5 ) )
97  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2)); } // Disc 3
98  // Disc 4 TEC2TEC
99  else if ( ( (theStep->GetPreStepPoint()->GetPosition().z() > 1682.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1802.5 )
100  || (theStep->GetPreStepPoint()->GetPosition().z() < -1682.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1802.5 ) )
101  && ( ( ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.285 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.295 )
102  || ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.84 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.86 )
103  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 3.63 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 3.66 )
104  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.20 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.23 )
105  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.76 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.80 ) ) ) )
106  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2*0.2)); } // Disc 4 TEC2TEC
107  // Disc 4
108  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 1682.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1802.5 )
109  || (theStep->GetPreStepPoint()->GetPosition().z() < -1682.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1802.5 ) )
110  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2)); } // Disc 4
111  // Disc 5 TEC2TEC
112  else if ( ( ( theStep->GetPreStepPoint()->GetPosition().z() > 1822.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1942.5 )
113  || ( theStep->GetPreStepPoint()->GetPosition().z() < -1822.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1942.5 ) )
114  && ( ( ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.285 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.295 )
115  || ( theStep->GetPreStepPoint()->GetPosition().phi() > 1.84 && theStep->GetPreStepPoint()->GetPosition().phi() < 1.86 )
116  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 3.63 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 3.66 )
117  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.20 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.23 )
118  || ( theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI > 5.76 && theStep->GetPreStepPoint()->GetPosition().phi() + 2.0*M_PI < 5.80 ) ) ) )
119  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2*0.2*0.2)); } // Disc 5 TEC2TEC
120  // Disc 5
121  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 1822.5 && theStep->GetPreStepPoint()->GetPosition().z() < 1942.5 )
122  || (theStep->GetPreStepPoint()->GetPosition().z() < -1822.5 && theStep->GetPreStepPoint()->GetPosition().z() > -1942.5 ) )
123  { theStep->AddTotalEnergyDeposit(EnergyLoss); } // Disc 5
124  // Disc 6
125  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 1997.5 && theStep->GetPreStepPoint()->GetPosition().z() < 2117.5 )
126  || (theStep->GetPreStepPoint()->GetPosition().z() < -1997.5 && theStep->GetPreStepPoint()->GetPosition().z() > -2117.5 ) )
127  { theStep->AddTotalEnergyDeposit(EnergyLoss); } // Disc 6
128  // Disc 7
129  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 2187.5 && theStep->GetPreStepPoint()->GetPosition().z() < 2307.5 )
130  || (theStep->GetPreStepPoint()->GetPosition().z() < -2187.5 && theStep->GetPreStepPoint()->GetPosition().z() > -2307.5 ) )
131  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2)); } // Disc 7
132  // Disc 8
133  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 2392.5 && theStep->GetPreStepPoint()->GetPosition().z() < 2512.5 )
134  || (theStep->GetPreStepPoint()->GetPosition().z() < -2392.5 && theStep->GetPreStepPoint()->GetPosition().z() > -2512.5 ) )
135  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2)); } // Disc 8
136  // Disc 9
137  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > 2607.5 && theStep->GetPreStepPoint()->GetPosition().z() < 2727.5 )
138  || (theStep->GetPreStepPoint()->GetPosition().z() < -2607.5 && theStep->GetPreStepPoint()->GetPosition().z() > -2727.5 ) )
139  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2*0.2)); } // Disc 9
140  // Beams in Barrel
141  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > -870.0 && theStep->GetPreStepPoint()->GetPosition().z() < 1050.0) &&
142  (theStep->GetPreStepPoint()->GetPosition().perp() > 500.0 && theStep->GetPreStepPoint()->GetPosition().perp() < 630.0) )
143  { theStep->AddTotalEnergyDeposit(EnergyLoss/(0.2*0.2)); } // Beams in the Barrel
144  else
145  {
146  // apparently we are not in a detector which should be hit by a LaserBeam
147  // therefore set the EnergyDeposit to zero and do not create a SimHit
148  theStep->ResetTotalEnergyDeposit();
149  }
150  }
151  }
152  // kill the photon if it goes through a module in the outer barrel detector. In practice on the back
153  // of a module is a thin layer of Aluminium that absorbs the photons, so hits will only be created in
154  // the first layer of the TOB. In the current geometry this Aluminium layer is not included. This should
155  // also avoid unwanted reflections (which then create hits in the TIB at the wrong positions)
156  else if ( ( (theStep->GetPreStepPoint()->GetMaterial()->GetName() == "TOB_Wafer") &&
157  (theStep->GetPostStepPoint()->GetMaterial()->GetName() == "Air") ) )
158  {
159  LogDebug("LaserAlignmentSimulationStepping") << " AC1CMS: stepping aborted! ";
160  theTrack->SetTrackStatus(fStopAndKill);
161  }
162  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > -870.0 &&
163  theStep->GetPreStepPoint()->GetPosition().z() < 1050.0) &&
164  (theStep->GetPreStepPoint()->GetPosition().perp() > 630.0) )
165  {
166  LogDebug("LaserAlignmentSimulationStepping") << " AC1CMS: stepping aborted! ";
167  theTrack->SetTrackStatus(fStopAndKill);
168  }
169  // do the same for photons that a) go through a module in the inner barrel detector or b) are reflected
170  // at the surface of a TIB module. The photons in case b) can create hits in the TOB at the wrong z
171  // positions :-(
172  else if ( ( (theStep->GetPreStepPoint()->GetMaterial()->GetName() == "TIB_Wafer") &&
173  (theStep->GetPostStepPoint()->GetMaterial()->GetName() == "Air") ) )
174  {
175  LogDebug("LaserAlignmentSimulationStepping") << " AC1CMS: stepping aborted! ";
176  theTrack->SetTrackStatus(fStopAndKill);
177  }
178  else if ( (theStep->GetPreStepPoint()->GetPosition().z() > -870.0 &&
179  theStep->GetPreStepPoint()->GetPosition().z() < 1050.0) &&
180  (theStep->GetPreStepPoint()->GetPosition().perp() < 500.0) )
181  {
182  LogDebug("LaserAlignmentSimulationStepping") << " AC1CMS: stepping aborted! ";
183  theTrack->SetTrackStatus(fStopAndKill);
184  }
185  // avoid reflections at Disc 1 of TEC- which enter again the Barrel Volume. These Photons
186  // create hits at the wrong positions in TIB and TOB
187  else if ( ( ( (theStep->GetPreStepPoint()->GetMaterial()->GetName() == "TEC_Wafer") &&
188  (theStep->GetPostStepPoint()->GetMaterial()->GetName() == "T_Air") ) ||
189  ( (theStep->GetPreStepPoint()->GetMaterial()->GetName() == "TEC_Wafer") &&
190  (theStep->GetPostStepPoint()->GetMaterial()->GetName() == "Air") ) ) &&
191  (theStep->GetPreStepPoint()->GetMomentum().z() != theStep->GetPostStepPoint()->GetMomentum().z() ) &&
192  (theStep->GetPostStepPoint()->GetPosition().z() == -1137.25) )
193  {
194  LogDebug("LaserAlignmentSimulationStepping") << " AC1CMS: stepping aborted! photon in wrong direction";
195  theTrack->SetTrackStatus(fStopAndKill);
196  }
197  // kill photons in the barrel which go in the wrong (i.e. +z) direction; they create unwanted hits
198  // due to reflections ...
199  else if ( ( theStep->GetPostStepPoint()->GetPosition().z() > -1100.0 )
200  && ( theStep->GetPostStepPoint()->GetPosition().z() < 1100.0 )
201  && ( theStep->GetPostStepPoint()->GetMomentumDirection().z() > 0.8 ) )
202  {
203  LogDebug("LaserAlignmentSimulationStepping") << " AC1CMS: stepping aborted! photon in wrong direction";
204  theTrack->SetTrackStatus(fStopAndKill);
205  }
206  else
207  {
208  LogDebug("LaserAlignmentSimulationStepping") << " AC1CMS: stepping continuous ... ";
209  }
210  // ***********************************************************************************************************
211 
212 
213  // check if it is alive
214  if ( theTrack->GetTrackStatus() != fAlive )
215  {
216  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: Track is not alive! -> return ";
217  return;
218  }
219 
220  // check if it is a primary
221  if ( theTrack->GetParentID() != 0 )
222  {
223  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: Track is not a primary! -> return ";
224  return;
225  }
226 
227  // check if it is a optical photon
228  if ( theDebugLevel >= 4 )
229  {
230  G4ParticleDefinition * theOpticalType = theTrack->GetDefinition();
231  if ( theOpticalType == G4OpticalPhoton::OpticalPhotonDefinition() )
232  {
233  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: Optical Photon found! ";
234  }
235 
236  // check in which volume it is
237  G4StepPoint * thePreStepPoint = theStep->GetPreStepPoint();
238  G4VPhysicalVolume * thePreStepPhysicalVolume = thePreStepPoint->GetPhysicalVolume();
239  G4String thePreStepPhysicalVolumeName = thePreStepPhysicalVolume->GetName();
240  G4Material * thePreStepMaterial = thePreStepPoint->GetMaterial();
241 
242  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: the PreStep Position = " << thePreStepPoint->GetPosition()
243  << "\n<LaserSteppingAction::UserSteppingAction(const G4Step *)>: the PreStep Physical Volume = " << thePreStepPhysicalVolumeName
244  << "\n<LaserSteppingAction::UserSteppingAction(const G4Step *)>: the PreStep Material =" << thePreStepMaterial->GetName();
245 
246  G4StepPoint * thePostStepPoint = theStep->GetPostStepPoint();
247  G4VPhysicalVolume * thePostStepPhysicalVolume = thePostStepPoint->GetPhysicalVolume();
248  G4String thePostStepPhysicalVolumeName = thePostStepPhysicalVolume->GetName();
249  G4Material * thePostStepMaterial = thePostStepPoint->GetMaterial();
250 
251  LogDebug("LaserAlignmentSimulationStepping") << "<LaserSteppingAction::UserSteppingAction(const G4Step *)>: the PostStep Position = " << thePostStepPoint->GetPosition()
252  << "\n<LaserSteppingAction::UserSteppingAction(const G4Step *)>: the PostStep Physical Volume = " << thePostStepPhysicalVolumeName
253  << "\n<LaserSteppingAction::UserSteppingAction(const G4Step *)>: the PostStep Material = " << thePostStepMaterial->GetName();
254  }
255 }
#define LogDebug(id)
virtual void UserSteppingAction(const G4Step *myStep)
stepping action: set energydeposit when a photon is absorbed in a Si module
virtual ~LaserSteppingAction()
destructor
#define M_PI
LaserSteppingAction(edm::ParameterSet const &theConf)
constructor