CMS 3D CMS Logo

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