CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/Alignment/LaserAlignmentSimulation/plugins/LaserAlignmentSimulation.cc

Go to the documentation of this file.
00001 
00009 #include "Alignment/LaserAlignmentSimulation/plugins/LaserAlignmentSimulation.h"
00010 #include "SimG4Core/Notification/interface/BeginOfRun.h" 
00011 #include "SimG4Core/Notification/interface/EndOfRun.h" 
00012 #include "SimG4Core/Notification/interface/BeginOfEvent.h" 
00013 #include "SimG4Core/Notification/interface/EndOfEvent.h" 
00014 #include "G4SDManager.hh" 
00015 #include "G4Step.hh" 
00016 #include "G4Timer.hh" 
00017 #include "G4VProcess.hh" 
00018 
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 
00021 
00022 #include "SimG4CMS/Tracker/interface/TkAccumulatingSensitiveDetector.h"
00023 #include "G4StepPoint.hh" 
00024 
00025 
00026 
00027 LaserAlignmentSimulation::LaserAlignmentSimulation(edm::ParameterSet const& theConf) 
00028   : theDebugLevel(theConf.getUntrackedParameter<int>("DebugLevel",0)),
00029     theEnergyLossScalingFactor(theConf.getUntrackedParameter<double>("EnergyLossScalingFactor",1.0)),
00030     theMPDebug(theConf.getUntrackedParameter<int>("MaterialPropertiesDebugLevel",0)),
00031     theSiAbsLengthScale(theConf.getUntrackedParameter<double>("SiAbsorptionLengthScalingFactor",1.0)),
00032     theTimer(), 
00033     theMaterialProperties(),
00034     thePrimaryGenerator(), theSteppingAction(),
00035     theBarrelHits(0), theEndcapHits(0),
00036                 theParameterSet(theConf)
00037 {
00038 
00039   // make some noise
00040   edm::LogInfo("SimLaserAlignmentSimulation") << " *****     AC1CMS: Configuration from ParameterSet      ***** " 
00041                                     << "\n  AC1CMS: theDebugLevel               = " << theDebugLevel 
00042                                     << "\n  AC1CMS: theEnergyLossScalingFactor  = " << theEnergyLossScalingFactor 
00043                                     << "\n  AC1CMS: theMPDebugLevel             = " << theMPDebug
00044                                     << "\n  AC1CMS: theSiAbsLengthScalingFactor = " << theSiAbsLengthScale;
00045 
00046   // declare timer
00047   theTimer = new G4Timer;
00048 }
00049 
00050 LaserAlignmentSimulation::~LaserAlignmentSimulation() 
00051 {
00052   if ( theMaterialProperties != 0 )        { delete theMaterialProperties; }
00053   if ( theSteppingAction != 0 )            { delete theSteppingAction; }
00054   if ( thePrimaryGenerator != 0 )          { delete thePrimaryGenerator; }
00055   if ( theTimer != 0 )                     { delete theTimer; }
00056 }
00057 
00058 void LaserAlignmentSimulation::update(const BeginOfRun * myRun)
00059 {
00060   LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfRun * myRun)>"
00061                                 << "\n *****     AC1CMS: Start of Run: " << (*myRun)()->GetRunID() << "     ***** ";
00062 
00063   // start timer
00064   theTimer->Start();
00065 
00066 
00067   // the PrimaryGeneratorAction: defines the used particlegun for the Laser events
00068   thePrimaryGenerator = new LaserPrimaryGeneratorAction(theParameterSet);
00069 
00070   // the UserSteppingAction: at the moment this prints only some information
00071   theSteppingAction = new LaserSteppingAction(theParameterSet);
00072 
00073   // construct your own material properties for setting refractionindex and so on
00074   theMaterialProperties = new MaterialProperties(theMPDebug, theSiAbsLengthScale);
00075 
00076   // list the tree of sensitive detectors
00077   if (theDebugLevel >= 1)
00078     {
00079       G4SDManager * theSDManager = G4SDManager::GetSDMpointer();
00080       theSDManager->ListTree();
00081     }
00082 }
00083 
00084 void LaserAlignmentSimulation::update(const BeginOfEvent * myEvent) 
00085 {
00086   LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const BeginOfEvent * myEvent)>"
00087                                 << "\n AC1CMS: Event number = " << (*myEvent)()->GetEventID();
00088 
00089   // some statistics for this event
00090   theBarrelHits = 0;
00091   theEndcapHits = 0;
00092 
00093   // generate the Primaries
00094   thePrimaryGenerator->GeneratePrimaries((G4Event*)(*myEvent)());
00095 }
00096 
00097 void LaserAlignmentSimulation::update(const BeginOfTrack * myTrack)
00098 {
00099 }
00100 
00101 void LaserAlignmentSimulation::update(const G4Step * myStep) 
00102 {
00103   LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step * myStep)>";
00104 
00105   G4Step * theStep = const_cast<G4Step*>(myStep);
00106 
00107   // do the LaserSteppingAction
00108   theSteppingAction->UserSteppingAction(theStep);
00109 
00110   // Trigger sensitive detector manually since photon is absorbed
00111   if ( ( theStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName()== "OpAbsorption" ) )
00112     {
00113       LogDebug("SimLaserAlignmentSimulationStepping") << "<LaserAlignmentSimulation::update(const G4Step*)>: Photon was absorbed! ";
00114       
00115       
00116       if ( theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector() )
00117         {
00118           LogDebug("SimLaserAlignmentSimulationStepping") << " AC1CMS: Setting the EnergyLoss to " << theStep->GetTotalEnergyDeposit() 
00119                                                 << "\n AC1CMS: The z position is " << theStep->GetPreStepPoint()->GetPosition().z()
00120                                                 << "\n AC1CMS: the Sensitive Detector: " 
00121                                                 << theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName()
00122                                                 << "\n AC1CMS: the Material: " << theStep->GetPreStepPoint()->GetMaterial()->GetName()
00123                                                 << "\n AC1CMS: the Logical Volume: " 
00124                                                 << theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
00125 
00126           if (theStep->GetTotalEnergyDeposit() > 0.0)
00127             {
00128               // process a hit
00129               TkAccumulatingSensitiveDetector * theSD = (TkAccumulatingSensitiveDetector*)
00130                 (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector());
00131               
00132               theSD->ProcessHits(theStep, ((G4TouchableHistory *)(theStep->GetPreStepPoint()->GetTouchable())));
00133 
00134 
00135               // some statistics for this event
00136               if ( ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TECModule3RphiActive" ) || 
00137                    ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TECModule5RphiActive" ) )
00138                 {
00139                   theEndcapHits++;
00140                 }
00141               else if ( ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TOBActiveSter0" ) ||
00142                         ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TOBActiveRphi0" ) ||
00143                         ( theStep->GetPostStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName() == "TIBActiveRphi2" ) )
00144                 {
00145                   theBarrelHits++;
00146                 }
00147             }
00148         }
00149       else 
00150         {
00151           LogDebug("SimLaserAlignmentSimulationStepping") << " AC1CMS: No SensitiveDetector available for this Step ... No Hit created :-( "
00152                                                 << "\n AC1CMS: The Material was: " << theStep->GetPreStepPoint()->GetMaterial()->GetName(); 
00153         }
00154     }
00155 }
00156 
00157 void LaserAlignmentSimulation::update(const EndOfTrack * myTrack)
00158 {
00159 }
00160 
00161 void LaserAlignmentSimulation::update(const EndOfEvent * myEvent) 
00162 {
00163   LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfEvent * myEvent)>"
00164                                 << "\n AC1CMS: End of Event " << (*myEvent)()->GetEventID();
00165 
00166   // some statistics for this event
00167   edm::LogInfo("SimLaserAlignmentSimulation") << " *** Number of Hits: " << theBarrelHits << " / " << theEndcapHits
00168                                     << " (Barrel / Endcaps) *** ";
00169 }
00170 
00171 void LaserAlignmentSimulation::update(const EndOfRun * myRun)
00172 {
00173   LogDebug("SimLaserAlignmentSimulation") << "<LaserAlignmentSimulation::update(const EndOfRun * myRun)>";
00174 
00175   // stop timer
00176   theTimer->Stop();
00177   edm::LogInfo("SimLaserAlignmentSimulation") << " AC1CMS: Number of Events = " << (*myRun)()->GetNumberOfEventToBeProcessed()
00178                                     << " " << *theTimer 
00179                                     << " *****     AC1CMS: End of Run: " << (*myRun)()->GetRunID() << "     ***** ";
00180 }
00181 
00182 // register a SimWatcher to get the Observer signals from OscarProducer
00183 #include "SimG4Core/Watcher/interface/SimWatcherFactory.h"
00184 
00185 
00186 DEFINE_SIMWATCHER (LaserAlignmentSimulation);