File indexing completed on 2024-04-06 11:56:27
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "Alignment/LaserAlignmentSimulation/plugins/LaserAlignmentSimulation.h"
0011 #include "G4SDManager.hh"
0012 #include "G4Step.hh"
0013 #include "G4Timer.hh"
0014 #include "G4VProcess.hh"
0015 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
0016 #include "SimG4Core/Notification/interface/BeginOfRun.h"
0017 #include "SimG4Core/Notification/interface/EndOfEvent.h"
0018 #include "SimG4Core/Notification/interface/EndOfRun.h"
0019
0020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0021
0022 #include "G4StepPoint.hh"
0023 #include "SimG4CMS/Tracker/interface/TkAccumulatingSensitiveDetector.h"
0024
0025 LaserAlignmentSimulation::LaserAlignmentSimulation(edm::ParameterSet const &theConf)
0026 : theDebugLevel(theConf.getUntrackedParameter<int>("DebugLevel", 0)),
0027 theEnergyLossScalingFactor(theConf.getUntrackedParameter<double>("EnergyLossScalingFactor", 1.0)),
0028 theMPDebug(theConf.getUntrackedParameter<int>("MaterialPropertiesDebugLevel", 0)),
0029 theSiAbsLengthScale(theConf.getUntrackedParameter<double>("SiAbsorptionLengthScalingFactor", 1.0)),
0030 theTimer(),
0031 theMaterialProperties(),
0032 thePrimaryGenerator(),
0033 theSteppingAction(),
0034 theBarrelHits(0),
0035 theEndcapHits(0),
0036 theParameterSet(theConf) {
0037
0038 edm::LogInfo("SimLaserAlignmentSimulation")
0039 << " ***** AC1CMS: Configuration from ParameterSet ***** "
0040 << "\n AC1CMS: theDebugLevel = " << theDebugLevel
0041 << "\n AC1CMS: theEnergyLossScalingFactor = " << theEnergyLossScalingFactor
0042 << "\n AC1CMS: theMPDebugLevel = " << theMPDebug
0043 << "\n AC1CMS: theSiAbsLengthScalingFactor = " << theSiAbsLengthScale;
0044
0045
0046 theTimer = new G4Timer;
0047 }
0048
0049 LaserAlignmentSimulation::~LaserAlignmentSimulation() {
0050 if (theMaterialProperties != nullptr) {
0051 delete theMaterialProperties;
0052 }
0053 if (theSteppingAction != nullptr) {
0054 delete theSteppingAction;
0055 }
0056 if (thePrimaryGenerator != nullptr) {
0057 delete thePrimaryGenerator;
0058 }
0059 if (theTimer != nullptr) {
0060 delete theTimer;
0061 }
0062 }
0063
0064 void LaserAlignmentSimulation::update(const BeginOfRun *myRun) {
0065 LogDebug("SimLaserAlignmentSimulation")
0066 << "<LaserAlignmentSimulation::update(const BeginOfRun * myRun)>"
0067 << "\n ***** AC1CMS: Start of Run: " << (*myRun)()->GetRunID() << " ***** ";
0068
0069
0070 theTimer->Start();
0071
0072
0073
0074 thePrimaryGenerator = new LaserPrimaryGeneratorAction(theParameterSet);
0075
0076
0077 theSteppingAction = new LaserSteppingAction(theParameterSet);
0078
0079
0080
0081 theMaterialProperties = new MaterialProperties(theMPDebug, theSiAbsLengthScale);
0082
0083
0084 if (theDebugLevel >= 1) {
0085 G4SDManager *theSDManager = G4SDManager::GetSDMpointer();
0086 theSDManager->ListTree();
0087 }
0088 }
0089
0090 void LaserAlignmentSimulation::update(const BeginOfEvent *myEvent) {
0091 LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfEvent * myEvent)>"
0092 << "\n AC1CMS: Event number = " << (*myEvent)()->GetEventID();
0093
0094
0095 theBarrelHits = 0;
0096 theEndcapHits = 0;
0097
0098
0099 thePrimaryGenerator->GeneratePrimaries((G4Event *)(*myEvent)());
0100 }
0101
0102 void LaserAlignmentSimulation::update(const BeginOfTrack *myTrack) {}
0103
0104 void LaserAlignmentSimulation::update(const G4Step *myStep) {
0105 LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step * myStep)>";
0106
0107 G4Step *theStep = const_cast<G4Step *>(myStep);
0108
0109
0110 theSteppingAction->UserSteppingAction(theStep);
0111
0112
0113 if ((theStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName() == "OpAbsorption")) {
0114 LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step*)>: Photon was "
0115 "absorbed! ";
0116
0117 if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()) {
0118 LogDebug("SimLaserAlignmentSimulationStepping")
0119 << " AC1CMS: Setting the EnergyLoss to " << theStep->GetTotalEnergyDeposit()
0120 << "\n AC1CMS: The z position is " << theStep->GetPreStepPoint()->GetPosition().z()
0121 << "\n AC1CMS: the Sensitive Detector: "
0122 << theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName()
0123 << "\n AC1CMS: the Material: " << theStep->GetPreStepPoint()->GetMaterial()->GetName()
0124 << "\n AC1CMS: the Logical Volume: "
0125 << theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
0126
0127 if (theStep->GetTotalEnergyDeposit() > 0.0) {
0128
0129 TkAccumulatingSensitiveDetector *theSD =
0130 (TkAccumulatingSensitiveDetector
0131 *)(theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector());
0132
0133 theSD->ProcessHits(theStep, ((G4TouchableHistory *)(theStep->GetPreStepPoint()->GetTouchable())));
0134
0135
0136 if ((theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
0137 "TECModule3RphiActive") ||
0138 (theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
0139 "TECModule5RphiActive")) {
0140 theEndcapHits++;
0141 } else if ((theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
0142 "TOBActiveSter0") ||
0143 (theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
0144 "TOBActiveRphi0") ||
0145 (theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() ==
0146 "TIBActiveRphi2")) {
0147 theBarrelHits++;
0148 }
0149 }
0150 } else {
0151 LogDebug("SimLaserAlignmentSimulationStepping")
0152 << " AC1CMS: No SensitiveDetector available for this Step ... No Hit "
0153 "created :-( "
0154 << "\n AC1CMS: The Material was: " << theStep->GetPreStepPoint()->GetMaterial()->GetName();
0155 }
0156 }
0157 }
0158
0159 void LaserAlignmentSimulation::update(const EndOfTrack *myTrack) {}
0160
0161 void LaserAlignmentSimulation::update(const EndOfEvent *myEvent) {
0162 LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfEvent * myEvent)>"
0163 << "\n AC1CMS: End of Event " << (*myEvent)()->GetEventID();
0164
0165
0166 edm::LogInfo("SimLaserAlignmentSimulation")
0167 << " *** Number of Hits: " << theBarrelHits << " / " << theEndcapHits << " (Barrel / Endcaps) *** ";
0168 }
0169
0170 void LaserAlignmentSimulation::update(const EndOfRun *myRun) {
0171 LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfRun * myRun)>";
0172
0173
0174 theTimer->Stop();
0175 edm::LogInfo("SimLaserAlignmentSimulation")
0176 << " AC1CMS: Number of Events = " << (*myRun)()->GetNumberOfEventToBeProcessed() << " " << *theTimer
0177 << " ***** AC1CMS: End of Run: " << (*myRun)()->GetRunID() << " ***** ";
0178 }
0179
0180
0181 #include "SimG4Core/Watcher/interface/SimWatcherFactory.h"
0182
0183 DEFINE_SIMWATCHER(LaserAlignmentSimulation);