00001 #include "G4Step.hh"
00002 #include "G4Electron.hh"
00003 #include "G4Positron.hh"
00004 #include "G4VProcess.hh"
00005 #include "G4SDManager.hh"
00006 #include "G4PrimaryParticle.hh"
00007 #include "G4PrimaryVertex.hh"
00008
00009 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00010 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
00011 #include "SimG4Core/Notification/interface/EndOfEvent.h"
00012
00013 #include "SimG4Core/GFlash/interface/GflashG4Watcher.h"
00014
00015 #include <TVector2.h>
00016
00017 const double radLength = 8.9;
00018 const double rMoliere = 21.9;
00019
00020
00021 GflashG4Watcher::GflashG4Watcher(const edm::ParameterSet& p) {
00022
00023 edm::ParameterSet myP = p.getParameter<edm::ParameterSet>("GflashG4Watcher");
00024 histFileName_ = myP.getParameter<std::string>("histFileName");
00025 recoEnergyScaleEB_ = myP.getParameter<double>("recoEnergyScaleEB");
00026 recoEnergyScaleEE_ = myP.getParameter<double>("recoEnergyScaleEE");
00027
00028
00029 histFile_ = new TFile(histFileName_.c_str(),"RECREATE");
00030
00031 TH1::AddDirectory(kTRUE);
00032
00033
00034 em_incE = new TH1F("em_incE","Incoming energy at Ecal;E (GeV);Number of Events",500,0.0,500.0);
00035 em_vtx_rho = new TH1F("em_vtx_rho","vertex position;#rho (cm);Number of Events",100,0.0,10.0);
00036 em_vtx_z = new TH1F("em_vtx_z","vertex position;z (cm);Number of Events",100,-10.0,10.0);
00037
00038 eb_ssp_rho = new TH1F("eb_ssp_rho","Shower starting position;#rho (cm);Number of Events",200,0.0,200.0);
00039 eb_hit_long = new TH1F("eb_hit_long","longitudinal hit position;shower depth (cm);Number of energy weighted hits",400,0.0,200.0);
00040 eb_hit_lat = new TH1F("eb_hit_lat","lateral hit position;arm (cm);Number of energy weighted hits",100,0.0,5.0);
00041 eb_hit_rz = new TH2F("eb_hit_rz","hit position along the shower direction;shower depth (cm);arm (cm)",400,0.0,200.0,100,0.0,5.0);
00042 eb_hit_long_sd = new TH1F("eb_hit_long_sd","longitudinal hit position in Sensitive Detector;shower depth (cm);Number of energy weighted hits",400,0.0,200.0);
00043 eb_hit_lat_sd = new TH1F("eb_hit_lat_sd","lateral hit position in Sensitive Detector;arm (cm);Number of energy weighted hits",100,0.0,5.0);
00044 eb_hit_rz_sd = new TH2F("eb_hit_rz_sd","hit position along the shower direction in Sensitive Detector;shower depth (cm);arm (cm)",400,0.0,200.0,100,0.0,5.0);
00045
00046 ee_ssp_z = new TH1F("ee_ssp_z","Shower starting position;z (cm);Number of Events",800,-400.0,400.0);
00047 ee_hit_long = new TH1F("ee_hit_long","longitudinal hit position;shower depth (cm);Number of energy weighted hits",800,0.0,400.0);
00048 ee_hit_lat = new TH1F("ee_hit_lat","lateral hit position;arm (cm);Number of energy weighted hits",100,0.0,5.0);
00049 ee_hit_rz = new TH2F("ee_hit_rz","hit position along the shower direction;shower depth (cm);arm (cm)",800,0.0,400.0,100,0.0,5.0);
00050 ee_hit_long_sd = new TH1F("ee_hit_long_sd","longitudinal hit position in Sensitive Detector;shower depth (cm);Number of energy weighted hits",800,0.0,400.0);
00051 ee_hit_lat_sd = new TH1F("ee_hit_lat_sd","lateral hit position in Sensitive Detector;arm (cm);Number of energy weighted hits",100,0.0,5.0);
00052 ee_hit_rz_sd = new TH2F("ee_hit_rz_sd","hit position along the shower direction in Sensitive Detector;shower depth (cm);arm (cm)",800,0.0,400.0,100,0.0,5.0);
00053
00054 }
00055
00056
00057 GflashG4Watcher::~GflashG4Watcher() {
00058 histFile_->cd();
00059 histFile_->Write();
00060 histFile_->Close();
00061 }
00062
00063
00064 void GflashG4Watcher::update(const BeginOfEvent* g4Event){
00065
00066 inc_flag = false;
00067
00068 const G4Event* evt = (*g4Event)();
00069 inc_vertex = evt->GetPrimaryVertex(0)->GetPosition();
00070 inc_position = inc_vertex;
00071 inc_direction = evt->GetPrimaryVertex(0)->GetPrimary(0)->GetMomentum().unit();
00072 inc_energy = evt->GetPrimaryVertex(0)->GetPrimary(0)->GetMomentum().mag();
00073 out_energy = 0;
00074
00075 em_incE->Fill(inc_energy/GeV);
00076 em_vtx_rho->Fill(inc_vertex.rho()/cm);
00077 em_vtx_z->Fill(inc_vertex.z()/cm);
00078
00079 if(std::abs(inc_direction.eta()) < 1.5) eb_ssp_rho->Fill(inc_position.rho()/cm);
00080 else ee_ssp_z->Fill(inc_position.z()/cm);
00081
00082 }
00083
00084
00085 void GflashG4Watcher::update(const EndOfEvent* g4Event){ }
00086
00087
00088 void GflashG4Watcher::update(const G4Step* aStep){
00089
00090 if(aStep == NULL) return;
00091
00092 double hitEnergy = aStep->GetTotalEnergyDeposit();
00093
00094 if(hitEnergy < 1.0e-6) return;
00095
00096 bool inEB = std::abs(inc_direction.eta()) < 1.5;
00097
00098 out_energy += hitEnergy;
00099
00100
00101 G4ThreeVector hitPosition = aStep->GetPreStepPoint()->GetPosition();
00102 G4ThreeVector diff = hitPosition - inc_position;
00103 double angle = diff.angle(inc_direction);
00104 double diff_z = std::abs(diff.mag() * std::cos(angle));
00105 double diff_r = std::abs(diff.mag() * std::sin(angle));
00106
00107 G4VSensitiveDetector* aSensitive = aStep->GetPreStepPoint()->GetSensitiveDetector();
00108
00109 if(inEB){
00110 hitEnergy *= recoEnergyScaleEB_;
00111 eb_hit_long->Fill(diff_z/cm,hitEnergy/GeV);
00112 eb_hit_lat->Fill(diff_r/cm,hitEnergy/GeV);
00113 eb_hit_rz->Fill(diff_z/cm,diff_r/cm,hitEnergy/GeV);
00114 if(aSensitive){
00115 eb_hit_long_sd->Fill(diff_z/cm,hitEnergy/GeV);
00116 eb_hit_lat_sd->Fill(diff_r/cm,hitEnergy/GeV);
00117 eb_hit_rz_sd->Fill(diff_z/cm,diff_r/cm,hitEnergy/GeV);
00118 }
00119 }
00120 else{
00121 hitEnergy *= recoEnergyScaleEE_;
00122 ee_hit_long->Fill(diff_z/cm,hitEnergy/GeV);
00123 ee_hit_lat->Fill(diff_r/cm,hitEnergy/GeV);
00124 ee_hit_rz->Fill(diff_z/cm,diff_r/cm,hitEnergy/GeV);
00125 if(aSensitive){
00126 ee_hit_long_sd->Fill(diff_z/cm,hitEnergy/GeV);
00127 ee_hit_lat_sd->Fill(diff_r/cm,hitEnergy/GeV);
00128 ee_hit_rz_sd->Fill(diff_z/cm,diff_r/cm,hitEnergy/GeV);
00129 }
00130 }
00131
00132
00133 }
00134
00135
00136
00137 #include "SimG4Core/Watcher/interface/SimWatcherFactory.h"
00138 #include "FWCore/PluginManager/interface/ModuleDef.h"
00139
00140
00141 DEFINE_SIMWATCHER(GflashG4Watcher);
00142