CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LaserAlignmentSimulation.cc
Go to the documentation of this file.
1 
14 #include "G4SDManager.hh"
15 #include "G4Step.hh"
16 #include "G4Timer.hh"
17 #include "G4VProcess.hh"
18 
20 
21 
23 #include "G4StepPoint.hh"
24 
25 
26 
28  : theDebugLevel(theConf.getUntrackedParameter<int>("DebugLevel",0)),
29  theEnergyLossScalingFactor(theConf.getUntrackedParameter<double>("EnergyLossScalingFactor",1.0)),
30  theMPDebug(theConf.getUntrackedParameter<int>("MaterialPropertiesDebugLevel",0)),
31  theSiAbsLengthScale(theConf.getUntrackedParameter<double>("SiAbsorptionLengthScalingFactor",1.0)),
32  theTimer(),
33  theMaterialProperties(),
34  thePrimaryGenerator(), theSteppingAction(),
35  theBarrelHits(0), theEndcapHits(0),
36  theParameterSet(theConf)
37 {
38 
39  // make some noise
40  edm::LogInfo("SimLaserAlignmentSimulation") << " ***** AC1CMS: Configuration from ParameterSet ***** "
41  << "\n AC1CMS: theDebugLevel = " << theDebugLevel
42  << "\n AC1CMS: theEnergyLossScalingFactor = " << theEnergyLossScalingFactor
43  << "\n AC1CMS: theMPDebugLevel = " << theMPDebug
44  << "\n AC1CMS: theSiAbsLengthScalingFactor = " << theSiAbsLengthScale;
45 
46  // declare timer
47  theTimer = new G4Timer;
48 }
49 
51 {
52  if ( theMaterialProperties != 0 ) { delete theMaterialProperties; }
53  if ( theSteppingAction != 0 ) { delete theSteppingAction; }
54  if ( thePrimaryGenerator != 0 ) { delete thePrimaryGenerator; }
55  if ( theTimer != 0 ) { delete theTimer; }
56 }
57 
59 {
60  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfRun * myRun)>"
61  << "\n ***** AC1CMS: Start of Run: " << (*myRun)()->GetRunID() << " ***** ";
62 
63  // start timer
64  theTimer->Start();
65 
66 
67  // the PrimaryGeneratorAction: defines the used particlegun for the Laser events
69 
70  // the UserSteppingAction: at the moment this prints only some information
72 
73  // construct your own material properties for setting refractionindex and so on
75 
76  // list the tree of sensitive detectors
77  if (theDebugLevel >= 1)
78  {
79  G4SDManager * theSDManager = G4SDManager::GetSDMpointer();
80  theSDManager->ListTree();
81  }
82 }
83 
85 {
86  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfEvent * myEvent)>"
87  << "\n AC1CMS: Event number = " << (*myEvent)()->GetEventID();
88 
89  // some statistics for this event
90  theBarrelHits = 0;
91  theEndcapHits = 0;
92 
93  // generate the Primaries
94  thePrimaryGenerator->GeneratePrimaries((G4Event*)(*myEvent)());
95 }
96 
98 {
99 }
100 
101 void LaserAlignmentSimulation::update(const G4Step * myStep)
102 {
103  LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step * myStep)>";
104 
105  G4Step * theStep = const_cast<G4Step*>(myStep);
106 
107  // do the LaserSteppingAction
109 
110  // Trigger sensitive detector manually since photon is absorbed
111  if ( ( theStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName()== "OpAbsorption" ) )
112  {
113  LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step*)>: Photon was absorbed! ";
114 
115 
116  if ( theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector() )
117  {
118  LogDebug("SimLaserAlignmentSimulationStepping") << " AC1CMS: Setting the EnergyLoss to " << theStep->GetTotalEnergyDeposit()
119  << "\n AC1CMS: The z position is " << theStep->GetPreStepPoint()->GetPosition().z()
120  << "\n AC1CMS: the Sensitive Detector: "
121  << theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName()
122  << "\n AC1CMS: the Material: " << theStep->GetPreStepPoint()->GetMaterial()->GetName()
123  << "\n AC1CMS: the Logical Volume: "
124  << theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
125 
126  if (theStep->GetTotalEnergyDeposit() > 0.0)
127  {
128  // process a hit
130  (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector());
131 
132  theSD->ProcessHits(theStep, ((G4TouchableHistory *)(theStep->GetPreStepPoint()->GetTouchable())));
133 
134 
135  // some statistics for this event
136  if ( ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TECModule3RphiActive" ) ||
137  ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TECModule5RphiActive" ) )
138  {
139  theEndcapHits++;
140  }
141  else if ( ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TOBActiveSter0" ) ||
142  ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TOBActiveRphi0" ) ||
143  ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TIBActiveRphi2" ) )
144  {
145  theBarrelHits++;
146  }
147  }
148  }
149  else
150  {
151  LogDebug("SimLaserAlignmentSimulationStepping") << " AC1CMS: No SensitiveDetector available for this Step ... No Hit created :-( "
152  << "\n AC1CMS: The Material was: " << theStep->GetPreStepPoint()->GetMaterial()->GetName();
153  }
154  }
155 }
156 
158 {
159 }
160 
162 {
163  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfEvent * myEvent)>"
164  << "\n AC1CMS: End of Event " << (*myEvent)()->GetEventID();
165 
166  // some statistics for this event
167  edm::LogInfo("SimLaserAlignmentSimulation") << " *** Number of Hits: " << theBarrelHits << " / " << theEndcapHits
168  << " (Barrel / Endcaps) *** ";
169 }
170 
172 {
173  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfRun * myRun)>";
174 
175  // stop timer
176  theTimer->Stop();
177  edm::LogInfo("SimLaserAlignmentSimulation") << " AC1CMS: Number of Events = " << (*myRun)()->GetNumberOfEventToBeProcessed()
178  << " " << *theTimer
179  << " ***** AC1CMS: End of Run: " << (*myRun)()->GetRunID() << " ***** ";
180 }
181 
182 // register a SimWatcher to get the Observer signals from OscarProducer
184 
185 
#define LogDebug(id)
#define DEFINE_SIMWATCHER(type)
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void UserSteppingAction(const G4Step *myStep)
stepping action: set energydeposit when a photon is absorbed in a Si module
virtual ~LaserAlignmentSimulation()
destructor
void update(const BeginOfRun *myRun)
observer for BeginOfRun
LaserAlignmentSimulation(edm::ParameterSet const &theConf)
constructor
LaserPrimaryGeneratorAction * thePrimaryGenerator
MaterialProperties * theMaterialProperties
LaserSteppingAction * theSteppingAction
void GeneratePrimaries(G4Event *myEvent)
call the corresponding GeneratePrimaries routines for both TEC&#39;s and the Barrel