CMS 3D CMS Logo

LaserAlignmentSimulation.cc
Go to the documentation of this file.
1 
11 #include "G4SDManager.hh"
12 #include "G4Step.hh"
13 #include "G4Timer.hh"
14 #include "G4VProcess.hh"
19 
21 
22 #include "G4StepPoint.hh"
24 
26  : theDebugLevel(theConf.getUntrackedParameter<int>("DebugLevel", 0)),
27  theEnergyLossScalingFactor(theConf.getUntrackedParameter<double>("EnergyLossScalingFactor", 1.0)),
28  theMPDebug(theConf.getUntrackedParameter<int>("MaterialPropertiesDebugLevel", 0)),
29  theSiAbsLengthScale(theConf.getUntrackedParameter<double>("SiAbsorptionLengthScalingFactor", 1.0)),
30  theTimer(),
31  theMaterialProperties(),
32  thePrimaryGenerator(),
33  theSteppingAction(),
34  theBarrelHits(0),
35  theEndcapHits(0),
36  theParameterSet(theConf) {
37  // make some noise
38  edm::LogInfo("SimLaserAlignmentSimulation")
39  << " ***** AC1CMS: Configuration from ParameterSet ***** "
40  << "\n AC1CMS: theDebugLevel = " << theDebugLevel
41  << "\n AC1CMS: theEnergyLossScalingFactor = " << theEnergyLossScalingFactor
42  << "\n AC1CMS: theMPDebugLevel = " << theMPDebug
43  << "\n AC1CMS: theSiAbsLengthScalingFactor = " << theSiAbsLengthScale;
44 
45  // declare timer
46  theTimer = new G4Timer;
47 }
48 
50  if (theMaterialProperties != nullptr) {
51  delete theMaterialProperties;
52  }
53  if (theSteppingAction != nullptr) {
54  delete theSteppingAction;
55  }
56  if (thePrimaryGenerator != nullptr) {
57  delete thePrimaryGenerator;
58  }
59  if (theTimer != nullptr) {
60  delete theTimer;
61  }
62 }
63 
65  LogDebug("SimLaserAlignmentSimulation")
66  << "<LaserAlignmentSimulation::update(const BeginOfRun * myRun)>"
67  << "\n ***** AC1CMS: Start of Run: " << (*myRun)()->GetRunID() << " ***** ";
68 
69  // start timer
70  theTimer->Start();
71 
72  // the PrimaryGeneratorAction: defines the used particlegun for the Laser
73  // events
75 
76  // the UserSteppingAction: at the moment this prints only some information
78 
79  // construct your own material properties for setting refractionindex and so
80  // on
82 
83  // list the tree of sensitive detectors
84  if (theDebugLevel >= 1) {
85  G4SDManager *theSDManager = G4SDManager::GetSDMpointer();
86  theSDManager->ListTree();
87  }
88 }
89 
91  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfEvent * myEvent)>"
92  << "\n AC1CMS: Event number = " << (*myEvent)()->GetEventID();
93 
94  // some statistics for this event
95  theBarrelHits = 0;
96  theEndcapHits = 0;
97 
98  // generate the Primaries
100 }
101 
103 
104 void LaserAlignmentSimulation::update(const G4Step *myStep) {
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  LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step*)>: Photon was "
115  "absorbed! ";
116 
117  if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()) {
118  LogDebug("SimLaserAlignmentSimulationStepping")
119  << " AC1CMS: Setting the EnergyLoss to " << theStep->GetTotalEnergyDeposit()
120  << "\n AC1CMS: The z position is " << theStep->GetPreStepPoint()->GetPosition().z()
121  << "\n AC1CMS: the Sensitive Detector: "
122  << theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName()
123  << "\n AC1CMS: the Material: " << theStep->GetPreStepPoint()->GetMaterial()->GetName()
124  << "\n AC1CMS: the Logical Volume: "
125  << theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
126 
127  if (theStep->GetTotalEnergyDeposit() > 0.0) {
128  // process a hit
131  *)(theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector());
132 
133  theSD->ProcessHits(theStep, ((G4TouchableHistory *)(theStep->GetPreStepPoint()->GetTouchable())));
134 
135  // some statistics for this event
136  if ((theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
137  "TECModule3RphiActive") ||
138  (theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
139  "TECModule5RphiActive")) {
140  theEndcapHits++;
141  } else if ((theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
142  "TOBActiveSter0") ||
143  (theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
144  "TOBActiveRphi0") ||
145  (theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
146  "TIBActiveRphi2")) {
147  theBarrelHits++;
148  }
149  }
150  } else {
151  LogDebug("SimLaserAlignmentSimulationStepping")
152  << " AC1CMS: No SensitiveDetector available for this Step ... No Hit "
153  "created :-( "
154  << "\n AC1CMS: The Material was: " << theStep->GetPreStepPoint()->GetMaterial()->GetName();
155  }
156  }
157 }
158 
160 
162  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfEvent * myEvent)>"
163  << "\n AC1CMS: End of Event " << (*myEvent)()->GetEventID();
164 
165  // some statistics for this event
166  edm::LogInfo("SimLaserAlignmentSimulation")
167  << " *** Number of Hits: " << theBarrelHits << " / " << theEndcapHits << " (Barrel / Endcaps) *** ";
168 }
169 
171  LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfRun * myRun)>";
172 
173  // stop timer
174  theTimer->Stop();
175  edm::LogInfo("SimLaserAlignmentSimulation")
176  << " AC1CMS: Number of Events = " << (*myRun)()->GetNumberOfEventToBeProcessed() << " " << *theTimer
177  << " ***** AC1CMS: End of Run: " << (*myRun)()->GetRunID() << " ***** ";
178 }
179 
180 // register a SimWatcher to get the Observer signals from OscarProducer
182 
LaserAlignmentSimulation::LaserAlignmentSimulation
LaserAlignmentSimulation(edm::ParameterSet const &theConf)
constructor
Definition: LaserAlignmentSimulation.cc:25
MessageLogger.h
LaserSteppingAction::UserSteppingAction
void UserSteppingAction(const G4Step *myStep) override
Definition: LaserSteppingAction.cc:20
MaterialProperties
Definition: MaterialProperties.h:15
LaserSteppingAction
Definition: LaserSteppingAction.h:16
TkAccumulatingSensitiveDetector
Definition: TkAccumulatingSensitiveDetector.h:31
LaserAlignmentSimulation::update
void update(const BeginOfRun *myRun) override
observer for BeginOfRun
Definition: LaserAlignmentSimulation.cc:64
LaserAlignmentSimulation::theDebugLevel
int theDebugLevel
Definition: LaserAlignmentSimulation.h:73
LaserAlignmentSimulation::theSiAbsLengthScale
double theSiAbsLengthScale
Definition: LaserAlignmentSimulation.h:76
EndOfEvent.h
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
LaserAlignmentSimulation
Definition: LaserAlignmentSimulation.h:42
EndOfTrack
Definition: EndOfTrack.h:6
TkAccumulatingSensitiveDetector.h
DEFINE_SIMWATCHER
#define DEFINE_SIMWATCHER(type)
Definition: SimWatcherFactory.h:13
BeginOfRun.h
TkAccumulatingSensitiveDetector::ProcessHits
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: TkAccumulatingSensitiveDetector.cc:104
BeginOfTrack
Definition: BeginOfTrack.h:6
LaserAlignmentSimulation::thePrimaryGenerator
LaserPrimaryGeneratorAction * thePrimaryGenerator
Definition: LaserAlignmentSimulation.h:81
EndOfEvent
Definition: EndOfEvent.h:6
ecalTB2006H4_GenSimDigiReco_cfg.myEvent
myEvent
Definition: ecalTB2006H4_GenSimDigiReco_cfg.py:7
LaserAlignmentSimulation.h
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
SimWatcherFactory.h
BeginOfEvent.h
createfilelist.int
int
Definition: createfilelist.py:10
LaserAlignmentSimulation::theBarrelHits
int theBarrelHits
Definition: LaserAlignmentSimulation.h:84
BeginOfEvent
Definition: BeginOfEvent.h:6
BeginOfRun
Definition: BeginOfRun.h:6
EndOfRun
Definition: EndOfRun.h:6
LaserAlignmentSimulation::~LaserAlignmentSimulation
~LaserAlignmentSimulation() override
destructor
Definition: LaserAlignmentSimulation.cc:49
EndOfRun.h
LaserAlignmentSimulation::theEnergyLossScalingFactor
double theEnergyLossScalingFactor
Definition: LaserAlignmentSimulation.h:74
LaserAlignmentSimulation::theMPDebug
int theMPDebug
Definition: LaserAlignmentSimulation.h:75
LaserAlignmentSimulation::theMaterialProperties
MaterialProperties * theMaterialProperties
Definition: LaserAlignmentSimulation.h:80
LaserAlignmentSimulation::theEndcapHits
int theEndcapHits
Definition: LaserAlignmentSimulation.h:85
LaserAlignmentSimulation::theParameterSet
edm::ParameterSet theParameterSet
Definition: LaserAlignmentSimulation.h:87
LaserPrimaryGeneratorAction::GeneratePrimaries
void GeneratePrimaries(G4Event *myEvent) override
Definition: LaserPrimaryGeneratorAction.cc:66
LaserAlignmentSimulation::theSteppingAction
LaserSteppingAction * theSteppingAction
Definition: LaserAlignmentSimulation.h:82
LaserAlignmentSimulation::theTimer
G4Timer * theTimer
Definition: LaserAlignmentSimulation.h:79
LaserPrimaryGeneratorAction
Definition: LaserPrimaryGeneratorAction.h:26