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 
22 
23 
25 #include "G4StepPoint.hh"
26 
27 
28 
30  : theDebugLevel(theConf.getUntrackedParameter<int>("DebugLevel",0)),
31  theEnergyLossScalingFactor(theConf.getUntrackedParameter<double>("EnergyLossScalingFactor",1.0)),
32  theMPDebug(theConf.getUntrackedParameter<int>("MaterialPropertiesDebugLevel",0)),
33  theSiAbsLengthScale(theConf.getUntrackedParameter<double>("SiAbsorptionLengthScalingFactor",1.0)),
34  theTimer(),
35  theMaterialProperties(),
36  thePrimaryGenerator(), theSteppingAction(),
37  theBarrelHits(0), theEndcapHits(0),
38  theParameterSet(theConf)
39 {
40 
41  // make some noise
42  edm::LogInfo("SimLaserAlignmentSimulation") << " ***** AC1CMS: Configuration from ParameterSet ***** "
43  << "\n AC1CMS: theDebugLevel = " << theDebugLevel
44  << "\n AC1CMS: theEnergyLossScalingFactor = " << theEnergyLossScalingFactor
45  << "\n AC1CMS: theMPDebugLevel = " << theMPDebug
46  << "\n AC1CMS: theSiAbsLengthScalingFactor = " << theSiAbsLengthScale;
47 
48  // declare timer
49  theTimer = new G4Timer;
50 }
51 
53 {
54  if ( theMaterialProperties != 0 ) { delete theMaterialProperties; }
55  if ( theSteppingAction != 0 ) { delete theSteppingAction; }
56  if ( thePrimaryGenerator != 0 ) { delete thePrimaryGenerator; }
57  if ( theTimer != 0 ) { delete theTimer; }
58 }
59 
61 {
62  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfRun * myRun)>"
63  << "\n ***** AC1CMS: Start of Run: " << (*myRun)()->GetRunID() << " ***** ";
64 
65  // start timer
66  theTimer->Start();
67 
68 
69  // the PrimaryGeneratorAction: defines the used particlegun for the Laser events
71 
72  // the UserSteppingAction: at the moment this prints only some information
74 
75  // construct your own material properties for setting refractionindex and so on
77 
78  // list the tree of sensitive detectors
79  if (theDebugLevel >= 1)
80  {
81  G4SDManager * theSDManager = G4SDManager::GetSDMpointer();
82  theSDManager->ListTree();
83  }
84 }
85 
87 {
88  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfEvent * myEvent)>"
89  << "\n AC1CMS: Event number = " << (*myEvent)()->GetEventID();
90 
91  // some statistics for this event
92  theBarrelHits = 0;
93  theEndcapHits = 0;
94 
95  // generate the Primaries
96  thePrimaryGenerator->GeneratePrimaries((G4Event*)(*myEvent)());
97 }
98 
100 {
101 }
102 
103 void LaserAlignmentSimulation::update(const G4Step * myStep)
104 {
105  LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step * myStep)>";
106 
107  G4Step * theStep = const_cast<G4Step*>(myStep);
108 
109  // do the LaserSteppingAction
111 
112  // Trigger sensitive detector manually since photon is absorbed
113  if ( ( theStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName()== "OpAbsorption" ) )
114  {
115  LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step*)>: Photon was absorbed! ";
116 
117 
118  if ( theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector() )
119  {
120  LogDebug("SimLaserAlignmentSimulationStepping") << " AC1CMS: Setting the EnergyLoss to " << theStep->GetTotalEnergyDeposit()
121  << "\n AC1CMS: The z position is " << theStep->GetPreStepPoint()->GetPosition().z()
122  << "\n AC1CMS: the Sensitive Detector: "
123  << theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName()
124  << "\n AC1CMS: the Material: " << theStep->GetPreStepPoint()->GetMaterial()->GetName()
125  << "\n AC1CMS: the Logical Volume: "
126  << theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
127 
128  if (theStep->GetTotalEnergyDeposit() > 0.0)
129  {
130  // process a hit
132  (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector());
133 
134  theSD->ProcessHits(theStep, ((G4TouchableHistory *)(theStep->GetPreStepPoint()->GetTouchable())));
135 
136 
137  // some statistics for this event
138  if ( ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TECModule3RphiActive" ) ||
139  ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TECModule5RphiActive" ) )
140  {
141  theEndcapHits++;
142  }
143  else if ( ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TOBActiveSter0" ) ||
144  ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TOBActiveRphi0" ) ||
145  ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TIBActiveRphi2" ) )
146  {
147  theBarrelHits++;
148  }
149  }
150  }
151  else
152  {
153  LogDebug("SimLaserAlignmentSimulationStepping") << " AC1CMS: No SensitiveDetector available for this Step ... No Hit created :-( "
154  << "\n AC1CMS: The Material was: " << theStep->GetPreStepPoint()->GetMaterial()->GetName();
155  }
156  }
157 }
158 
160 {
161 }
162 
164 {
165  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfEvent * myEvent)>"
166  << "\n AC1CMS: End of Event " << (*myEvent)()->GetEventID();
167 
168  // some statistics for this event
169  edm::LogInfo("SimLaserAlignmentSimulation") << " *** Number of Hits: " << theBarrelHits << " / " << theEndcapHits
170  << " (Barrel / Endcaps) *** ";
171 }
172 
174 {
175  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfRun * myRun)>";
176 
177  // stop timer
178  theTimer->Stop();
179  edm::LogInfo("SimLaserAlignmentSimulation") << " AC1CMS: Number of Events = " << (*myRun)()->GetNumberOfEventToBeProcessed()
180  << " " << *theTimer
181  << " ***** AC1CMS: End of Run: " << (*myRun)()->GetRunID() << " ***** ";
182 }
183 
184 // register a SimWatcher to get the Observer signals from OscarProducer
186 
187 
#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