Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
00013 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00014 #include "DataFormats/Math/interface/Point3D.h"
00015
00016 #include "SimG4CMS/Forward/interface/DoCastorAnalysis.h"
00017 #include "SimG4CMS/Forward/interface/CastorNumberingScheme.h"
00018
00019 #include "TFile.h"
00020 #include <cmath>
00021 #include <iostream>
00022 #include <iomanip>
00023 #include <memory>
00024 #include <vector>
00025
00026 #define debug 0
00027
00028 DoCastorAnalysis::DoCastorAnalysis(const edm::ParameterSet &p) {
00029
00030 edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("DoCastorAnalysis");
00031 verbosity = m_Anal.getParameter<int>("Verbosity");
00032
00033 TreeFileName = m_Anal.getParameter<std::string>("CastorTreeFileName");
00034
00035 if (verbosity > 0) {
00036
00037 std::cout<<std::endl;
00038 std::cout<<"============================================================================"<<std::endl;
00039 std::cout << "DoCastorAnalysis:: Initialized as observer"<< std::endl;
00040
00041 std::cout <<" Castor Tree will be created"<< std::endl;
00042 std::cout <<" Castor Tree will be in file: "<<TreeFileName<<std::endl;
00043 if(debug) getchar();
00044
00045 std::cout<<"============================================================================"<<std::endl;
00046 std::cout<<std::endl;
00047 }
00048
00049 std::cout << "DoCastorAnalysis: output event root file created"<< std::endl;
00050 TString treefilename = TreeFileName;
00051 CastorOutputEventFile = new TFile(treefilename,"RECREATE");
00052
00053 CastorTree = new TTree("Sim","Sim");
00054
00055 CastorTree->Branch("simhit_x","std::vector<double>",&psimhit_x);
00056 CastorTree->Branch("simhit_y","std::vector<double>",&psimhit_y);
00057 CastorTree->Branch("simhit_z","std::vector<double>",&psimhit_z);
00058
00059 CastorTree->Branch("simhit_eta","std::vector<double>",&psimhit_eta);
00060 CastorTree->Branch("simhit_phi","std::vector<double>",&psimhit_phi);
00061 CastorTree->Branch("simhit_energy","std::vector<double>",&psimhit_energy);
00062
00063
00064 CastorTree->Branch("simhit_sector","std::vector<int>",&psimhit_sector);
00065 CastorTree->Branch("simhit_module","std::vector<int>",&psimhit_module);
00066
00067 CastorTree->Branch("simhit_etot",&simhit_etot,"simhit_etot/D");
00068 }
00069
00070 DoCastorAnalysis::~DoCastorAnalysis() {
00071
00072
00073
00074 CastorOutputEventFile->cd();
00075
00076 CastorTree->Write("",TObject::kOverwrite);
00077 std::cout << "DoCastorAnalysis: Ntuple event written" << std::endl;
00078 if(debug) getchar();
00079 CastorOutputEventFile->Close();
00080 std::cout << "DoCastorAnalysis: Event file closed" << std::endl;
00081 if(debug) getchar();
00082
00083 if (verbosity > 0) {
00084 std::cout<<std::endl<<"DoCastorAnalysis: end of process"<<std::endl;
00085 if(debug) getchar();
00086 }
00087
00088 }
00089
00090
00091
00092 void DoCastorAnalysis::update(const BeginOfJob * job) {
00093
00094 std::cout << " Starting new job " << std::endl;
00095 }
00096
00097
00098
00099 void DoCastorAnalysis::update(const BeginOfRun * run) {
00100
00101 std::cout << std::endl << "DoCastorAnalysis: Starting Run"<< std::endl;
00102
00103
00104
00105
00106
00107 eventIndex = 1;
00108 }
00109
00110 void DoCastorAnalysis::update(const BeginOfEvent * evt) {
00111 std::cout << "DoCastorAnalysis: Processing Event Number: "<<eventIndex<< std::endl;
00112 eventIndex++;
00113 }
00114
00115
00116
00117
00118 void DoCastorAnalysis::update(const EndOfEvent * evt) {
00119
00120
00121
00122
00123 G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
00124
00125 int CAFIid = G4SDManager::GetSDMpointer()->GetCollectionID("CastorFI");
00126 CaloG4HitCollection* theCAFI = (CaloG4HitCollection*) allHC->GetHC(CAFIid);
00127
00128 CastorNumberingScheme *theCastorNumScheme = new CastorNumberingScheme();
00129
00130 unsigned int volumeID=0;
00131
00132
00133 int nentries = theCAFI->entries();
00134 if(debug) std::cout<<"nentries in CAFI: "<<nentries<<std::endl;
00135 if(debug) getchar();
00136
00137 psimhit_x=&simhit_x;
00138 psimhit_x->clear();
00139 psimhit_x->reserve(nentries);
00140
00141 psimhit_y=&simhit_y;
00142 psimhit_y->clear();
00143 psimhit_y->reserve(nentries);
00144
00145 psimhit_z=&simhit_z;
00146 psimhit_z->clear();
00147 psimhit_z->reserve(nentries);
00148
00149 psimhit_eta=&simhit_eta;
00150 psimhit_eta->clear();
00151 psimhit_eta->reserve(nentries);
00152
00153 psimhit_phi=&simhit_phi;
00154 psimhit_phi->clear();
00155 psimhit_phi->reserve(nentries);
00156
00157 psimhit_energy=&simhit_energy;
00158 psimhit_energy->clear();
00159 psimhit_energy->reserve(nentries);
00160
00161
00162
00163
00164
00165 psimhit_sector=&simhit_sector;
00166 psimhit_sector->clear();
00167 psimhit_sector->reserve(nentries);
00168
00169 psimhit_module=&simhit_module;
00170 psimhit_module->clear();
00171 psimhit_module->reserve(nentries);
00172
00173 simhit_etot = 0;
00174
00175 if (nentries > 0) {
00176
00177 for (int ihit = 0; ihit < nentries; ihit++) {
00178 CaloG4Hit* aHit = (*theCAFI)[ihit];
00179 volumeID = aHit->getUnitID();
00180
00181
00182 int zside,sector,zmodule;
00183
00184 theCastorNumScheme->unpackIndex(volumeID,zside,sector,zmodule);
00185
00186 double energy = aHit->getEnergyDeposit()/GeV;
00187
00188
00189 math::XYZPoint pos = aHit->getPosition();
00190 double theta = pos.theta();
00191 double eta = -log(tan(theta/2.));
00192 double phi = pos.phi();
00193
00194 psimhit_x->push_back(pos.x());
00195 psimhit_y->push_back(pos.y());
00196 psimhit_z->push_back(pos.z());
00197
00198 psimhit_eta->push_back(eta);
00199 psimhit_phi->push_back(phi);
00200 psimhit_energy->push_back(energy);
00201
00202
00203 psimhit_sector->push_back(sector);
00204 psimhit_module->push_back(zmodule);
00205
00206 simhit_etot+=energy;
00207
00208 if(debug) std::cout<<"hit "<<ihit+1<<" : x = "<<(*psimhit_x)[ihit]<<" , eta = "<<(*psimhit_eta)[ihit]
00209 <<" , phi = "<<(*psimhit_phi)[ihit]<<" , energy = "<<(*psimhit_energy)[ihit]<<std::endl;
00210 }
00211
00212
00213 if(debug) getchar();
00214 CastorTree->Fill();
00215
00216 }
00217 }
00218
00219 void DoCastorAnalysis::update(const EndOfRun * run) {;}
00220
00221 void DoCastorAnalysis::update(const G4Step * aStep) {;}
00222
00223